xref: /petsc/src/sys/utils/sortso.c (revision 3ca90d2d9fe4d5ec7086bd4aee14f89370d16392)
1676f2a66SJacob Faibussowitsch #include <petscsys.h>                /*I  "petscsys.h"  I*/
2676f2a66SJacob Faibussowitsch #include <petsc/private/petscimpl.h>
3676f2a66SJacob Faibussowitsch 
44d3610e3SJacob Faibussowitsch PETSC_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 
104d3610e3SJacob Faibussowitsch PETSC_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 
164d3610e3SJacob Faibussowitsch PETSC_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__)
28676f2a66SJacob Faibussowitsch PETSC_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 
366e260675SJacob Faibussowitsch PETSC_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 
47676f2a66SJacob Faibussowitsch PETSC_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 
53676f2a66SJacob Faibussowitsch PETSC_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 
60676f2a66SJacob Faibussowitsch PETSC_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 
66676f2a66SJacob Faibussowitsch PETSC_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
73676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void COPYSWAPPY(char *a, char *b, char *t, size_t size)
74676f2a66SJacob Faibussowitsch {
75676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
76676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
77676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(t, b, size);CHKERRABORT(PETSC_COMM_SELF,ierr);
78676f2a66SJacob Faibussowitsch   ierr = PetscMemmove(b, a, size);CHKERRABORT(PETSC_COMM_SELF,ierr);
79676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(a, t, size);CHKERRABORT(PETSC_COMM_SELF,ierr);
80676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
81676f2a66SJacob Faibussowitsch }
82676f2a66SJacob Faibussowitsch 
836e260675SJacob Faibussowitsch PETSC_STATIC_INLINE void COPYSWAPPY2(char *al, char *ar, size_t asize, char *bl, char *br, size_t bsize, char *t)
84676f2a66SJacob Faibussowitsch {
85676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
86676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
876e260675SJacob Faibussowitsch   ierr = PetscMemcpy(t, ar, asize);CHKERRABORT(PETSC_COMM_SELF,ierr);
886e260675SJacob Faibussowitsch   ierr = PetscMemmove(ar, al, asize);CHKERRABORT(PETSC_COMM_SELF,ierr);
896e260675SJacob Faibussowitsch   ierr = PetscMemcpy(al, t, asize);CHKERRABORT(PETSC_COMM_SELF,ierr);
906e260675SJacob Faibussowitsch   ierr = PetscMemcpy(t, br, bsize);CHKERRABORT(PETSC_COMM_SELF,ierr);
916e260675SJacob Faibussowitsch   ierr = PetscMemmove(br, bl, bsize);CHKERRABORT(PETSC_COMM_SELF,ierr);
926e260675SJacob Faibussowitsch   ierr = PetscMemcpy(bl, t, bsize);CHKERRABORT(PETSC_COMM_SELF,ierr);
93676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
94676f2a66SJacob Faibussowitsch }
95676f2a66SJacob Faibussowitsch 
96676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memcpy(char *dest, const char *src, size_t size)
97676f2a66SJacob Faibussowitsch {
98676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
99676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
100676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(dest, src, size);CHKERRABORT(PETSC_COMM_SELF,ierr);
101676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
102676f2a66SJacob Faibussowitsch }
103676f2a66SJacob Faibussowitsch 
104676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memcpy2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize)
105676f2a66SJacob Faibussowitsch {
106676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
107676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
108676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(adest, asrc, asize);CHKERRABORT(PETSC_COMM_SELF,ierr);
1096e260675SJacob Faibussowitsch   ierr = PetscMemcpy(bdest, bsrc, bsize);CHKERRABORT(PETSC_COMM_SELF,ierr);
110676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
111676f2a66SJacob Faibussowitsch }
112676f2a66SJacob Faibussowitsch 
113676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memmove(char *dest, const char *src, size_t size)
114676f2a66SJacob Faibussowitsch {
115676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
116676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
117676f2a66SJacob Faibussowitsch   ierr = PetscMemmove(dest, src, size);CHKERRABORT(PETSC_COMM_SELF,ierr);
118676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
119676f2a66SJacob Faibussowitsch }
120676f2a66SJacob Faibussowitsch 
121676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memmove2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize)
122676f2a66SJacob Faibussowitsch {
123676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
124676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
125676f2a66SJacob Faibussowitsch   ierr = PetscMemmove(adest, asrc, asize);CHKERRABORT(PETSC_COMM_SELF,ierr);
126676f2a66SJacob Faibussowitsch   ierr = PetscMemmove(bdest, bsrc, bsize);CHKERRABORT(PETSC_COMM_SELF,ierr);
127676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
128676f2a66SJacob Faibussowitsch }
129676f2a66SJacob Faibussowitsch #endif
130676f2a66SJacob Faibussowitsch 
131676f2a66SJacob 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] >
132676f2a66SJacob Faibussowitsch  x. Output also inclusive.
133676f2a66SJacob Faibussowitsch 
134676f2a66SJacob Faibussowitsch  NOTE: Both gallopsearch functions CANNOT distinguish between inserting AFTER the entire array vs inserting at entry n!! For example for an array:
135676f2a66SJacob Faibussowitsch 
136676f2a66SJacob Faibussowitsch    {0,1,2,3,4,5}
137676f2a66SJacob Faibussowitsch 
138676f2a66SJacob Faibussowitsch    when looking to insert "5" this routine will return *m = 6, but when looking to insert "6" it will ALSO return *m = 6.
139676f2a66SJacob Faibussowitsch   */
1404d3610e3SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscGallopSearchLeft_Private(const char *arr, size_t size, CompFunc cmp, void *ctx, PetscInt l, PetscInt r, const char *x, PetscInt *m)
141676f2a66SJacob Faibussowitsch {
142676f2a66SJacob Faibussowitsch   PetscInt last = l, k = 1, mid, cur = l+1;
143676f2a66SJacob Faibussowitsch 
144676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
145676f2a66SJacob Faibussowitsch   *m = l;
146*3ca90d2dSJacob Faibussowitsch   if (PetscUnlikelyDebug(r < l)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"r %" PetscInt_FMT " < l %" PetscInt_FMT " in PetscGallopSearchLeft",r,l);
1474d3610e3SJacob Faibussowitsch   if ((*cmp)(x, arr+r*size, ctx) >= 0) {*m = r; PetscFunctionReturn(0);}
1484d3610e3SJacob Faibussowitsch   if ((*cmp)(x, (arr)+l*size, ctx) < 0 || PetscUnlikely(!(r-l))) PetscFunctionReturn(0);
149676f2a66SJacob Faibussowitsch   while (PETSC_TRUE) {
150676f2a66SJacob Faibussowitsch     if (PetscUnlikely(cur > r)) {cur = r; break;}
1514d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr)+cur*size, ctx) < 0) break;
152676f2a66SJacob Faibussowitsch     last = cur;
153676f2a66SJacob Faibussowitsch     cur += (k <<= 1) + 1;
154676f2a66SJacob Faibussowitsch     ++k;
155676f2a66SJacob Faibussowitsch   }
156676f2a66SJacob Faibussowitsch   /* standard binary search but take last 0 mid 0 cur 1 into account*/
157676f2a66SJacob Faibussowitsch   while (cur > last + 1) {
158676f2a66SJacob Faibussowitsch     mid = last + ((cur - last) >> 1);
1594d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr)+mid*size, ctx) < 0) {
160676f2a66SJacob Faibussowitsch       cur = mid;
161676f2a66SJacob Faibussowitsch     } else {
162676f2a66SJacob Faibussowitsch       last = mid;
163676f2a66SJacob Faibussowitsch     }
164676f2a66SJacob Faibussowitsch   }
165676f2a66SJacob Faibussowitsch   *m = cur;
166676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
167676f2a66SJacob Faibussowitsch }
168676f2a66SJacob Faibussowitsch 
169676f2a66SJacob Faibussowitsch /* Start right look left. Looking for e.g. A[-1] in B or mergehi. l inclusive, r inclusive. Returns last m such that arr[m]
170676f2a66SJacob Faibussowitsch  < x. Output also inclusive */
1714d3610e3SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscGallopSearchRight_Private(const char *arr, size_t size, CompFunc cmp, void *ctx, PetscInt l, PetscInt r, const char *x, PetscInt *m)
172676f2a66SJacob Faibussowitsch {
173676f2a66SJacob Faibussowitsch   PetscInt last = r, k = 1, mid, cur = r-1;
174676f2a66SJacob Faibussowitsch 
175676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
176676f2a66SJacob Faibussowitsch   *m = r;
177*3ca90d2dSJacob Faibussowitsch   if (PetscUnlikelyDebug(r < l)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"r %" PetscInt_FMT " < l %" PetscInt_FMT " in PetscGallopSearchRight",r,l);
1784d3610e3SJacob Faibussowitsch   if ((*cmp)(x, arr+l*size, ctx) <= 0) {*m = l; PetscFunctionReturn(0);}
1794d3610e3SJacob Faibussowitsch   if ((*cmp)(x, (arr)+r*size, ctx) > 0 || PetscUnlikely(!(r-l))) PetscFunctionReturn(0);
180676f2a66SJacob Faibussowitsch   while (PETSC_TRUE) {
181676f2a66SJacob Faibussowitsch     if (PetscUnlikely(cur < l)) {cur = l; break;}
1824d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr)+cur*size, ctx) > 0) break;
183676f2a66SJacob Faibussowitsch     last = cur;
184676f2a66SJacob Faibussowitsch     cur -= (k <<= 1) + 1;
185676f2a66SJacob Faibussowitsch     ++k;
186676f2a66SJacob Faibussowitsch   }
187676f2a66SJacob Faibussowitsch   /* standard binary search but take last r-1 mid r-1 cur r-2 into account*/
188676f2a66SJacob Faibussowitsch   while (last > cur + 1) {
189676f2a66SJacob Faibussowitsch     mid = last - ((last - cur) >> 1);
1904d3610e3SJacob Faibussowitsch     if ((*cmp)(x, (arr)+mid*size, ctx) > 0) {
191676f2a66SJacob Faibussowitsch       cur = mid;
192676f2a66SJacob Faibussowitsch     } else {
193676f2a66SJacob Faibussowitsch       last = mid;
194676f2a66SJacob Faibussowitsch     }
195676f2a66SJacob Faibussowitsch   }
196676f2a66SJacob Faibussowitsch   *m = cur;
197676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
198676f2a66SJacob Faibussowitsch }
199676f2a66SJacob Faibussowitsch 
200676f2a66SJacob Faibussowitsch /* Mergesort where size of left half <= size of right half, so mergesort is done left to right. Arr should be pointer to
201676f2a66SJacob Faibussowitsch  complete array, left is first index of left array, mid is first index of right array, right is last index of right
202676f2a66SJacob Faibussowitsch  array */
2034d3610e3SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortMergeLo_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right)
204676f2a66SJacob Faibussowitsch {
205676f2a66SJacob Faibussowitsch   PetscInt       i = 0, j = mid, k = left, gallopleft = 0, gallopright = 0;
206676f2a66SJacob Faibussowitsch   const PetscInt llen = mid-left;
207676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
208676f2a66SJacob Faibussowitsch 
209676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
210676f2a66SJacob Faibussowitsch   Petsc_memcpy(tarr, arr+(left*size), llen*size);
211676f2a66SJacob Faibussowitsch   while ((i < llen) && (j <= right)) {
2124d3610e3SJacob Faibussowitsch     if ((*cmp)(tarr+(i*size), arr+(j*size), ctx) < 0) {
213676f2a66SJacob Faibussowitsch       Petsc_memcpy(arr+(k*size), tarr+(i*size), size);
214676f2a66SJacob Faibussowitsch       ++k;
215676f2a66SJacob Faibussowitsch       gallopright = 0;
216676f2a66SJacob Faibussowitsch       if (++i < llen && ++gallopleft >= MIN_GALLOP_GLOBAL) {
217676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
218676f2a66SJacob Faibussowitsch         do {
219676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
220676f2a66SJacob Faibussowitsch           /* search temp for right[j], can move up to that of temp into arr immediately */
2214d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(tarr, size, cmp, ctx, i, llen-1, arr+(j*size), &l1);CHKERRQ(ierr);
222676f2a66SJacob Faibussowitsch           diff1 = l1-i;
223676f2a66SJacob Faibussowitsch           Petsc_memcpy(arr+(k*size), tarr+(i*size), diff1*size);
224676f2a66SJacob Faibussowitsch           k += diff1;
225676f2a66SJacob Faibussowitsch           i = l1;
226676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
2274d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, size, cmp, ctx, j, right, tarr+(i*size), &l2);CHKERRQ(ierr);
228676f2a66SJacob Faibussowitsch           diff2 = l2-j;
229676f2a66SJacob Faibussowitsch           Petsc_memmove((arr)+k*size, (arr)+j*size, diff2*size);
230676f2a66SJacob Faibussowitsch           k += diff2;
231676f2a66SJacob Faibussowitsch           j = l2;
232676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
233676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
234676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
235676f2a66SJacob Faibussowitsch       }
236676f2a66SJacob Faibussowitsch     } else {
237676f2a66SJacob Faibussowitsch       Petsc_memmove(arr+(k*size), arr+(j*size), size);
238676f2a66SJacob Faibussowitsch       ++k;
239676f2a66SJacob Faibussowitsch       gallopleft = 0;
240676f2a66SJacob Faibussowitsch       if (++j <= right && ++gallopright >= MIN_GALLOP_GLOBAL) {
241676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
242676f2a66SJacob Faibussowitsch         do {
243676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
244676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
2454d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, size, cmp, ctx, j, right, tarr+(i*size), &l2);CHKERRQ(ierr);
246676f2a66SJacob Faibussowitsch           diff2 = l2-j;
247676f2a66SJacob Faibussowitsch           Petsc_memmove(arr+(k*size), arr+(j*size), diff2*size);
248676f2a66SJacob Faibussowitsch           k += diff2;
249676f2a66SJacob Faibussowitsch           j = l2;
250676f2a66SJacob Faibussowitsch           /* search temp for right[j], can copy up to that of temp into arr immediately */
2514d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(tarr, size, cmp, ctx, i, llen-1, arr+(j*size), &l1);CHKERRQ(ierr);
252676f2a66SJacob Faibussowitsch           diff1 = l1-i;
253676f2a66SJacob Faibussowitsch           Petsc_memcpy(arr+(k*size), tarr+(i*size), diff1*size);
254676f2a66SJacob Faibussowitsch           k += diff1;
255676f2a66SJacob Faibussowitsch           i = l1;
256676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
257676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
258676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
259676f2a66SJacob Faibussowitsch       }
260676f2a66SJacob Faibussowitsch     }
261676f2a66SJacob Faibussowitsch   }
262676f2a66SJacob Faibussowitsch   if (i<llen) {Petsc_memcpy(arr+(k*size), tarr+(i*size), (llen-i)*size);}
263676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
264676f2a66SJacob Faibussowitsch }
265676f2a66SJacob Faibussowitsch 
2664d3610e3SJacob Faibussowitsch PETSC_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)
267676f2a66SJacob Faibussowitsch {
268676f2a66SJacob Faibussowitsch   PetscInt       i = 0, j = mid, k = left, gallopleft = 0, gallopright = 0;
269676f2a66SJacob Faibussowitsch   const PetscInt llen = mid-left;
270676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
271676f2a66SJacob Faibussowitsch 
272676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
273676f2a66SJacob Faibussowitsch   Petsc_memcpy2(atarr, arr+(left*asize), llen*asize, btarr, barr+(left*bsize), llen*bsize);
274676f2a66SJacob Faibussowitsch   while ((i < llen) && (j <= right)) {
2754d3610e3SJacob Faibussowitsch     if ((*cmp)(atarr+(i*asize), arr+(j*asize), ctx) < 0) {
276676f2a66SJacob Faibussowitsch       Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), asize, barr+(k*bsize), btarr+(i*bsize), bsize);
277676f2a66SJacob Faibussowitsch       ++k;
278676f2a66SJacob Faibussowitsch       gallopright = 0;
279676f2a66SJacob Faibussowitsch       if (++i < llen && ++gallopleft >= MIN_GALLOP_GLOBAL) {
280676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
281676f2a66SJacob Faibussowitsch         do {
282676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
283676f2a66SJacob Faibussowitsch           /* search temp for right[j], can move up to that of temp into arr immediately */
2844d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(atarr, asize, cmp, ctx, i, llen-1, arr+(j*asize), &l1);CHKERRQ(ierr);
285676f2a66SJacob Faibussowitsch           diff1 = l1-i;
286676f2a66SJacob Faibussowitsch           Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), diff1*asize, barr+(k*bsize), btarr+(i*bsize), diff1*bsize);
287676f2a66SJacob Faibussowitsch           k += diff1;
288676f2a66SJacob Faibussowitsch           i = l1;
289676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
2904d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, j, right, atarr+(i*asize), &l2);CHKERRQ(ierr);
291676f2a66SJacob Faibussowitsch           diff2 = l2-j;
292676f2a66SJacob Faibussowitsch           Petsc_memmove2(arr+(k*asize), arr+(j*asize), diff2*asize, barr+(k*bsize), barr+(j*bsize), diff2*bsize);
293676f2a66SJacob Faibussowitsch           k += diff2;
294676f2a66SJacob Faibussowitsch           j = l2;
295676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
296676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
297676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
298676f2a66SJacob Faibussowitsch       }
299676f2a66SJacob Faibussowitsch     } else {
300676f2a66SJacob Faibussowitsch       Petsc_memmove2(arr+(k*asize), arr+(j*asize), asize, barr+(k*bsize), barr+(j*bsize), bsize);
301676f2a66SJacob Faibussowitsch       ++k;
302676f2a66SJacob Faibussowitsch       gallopleft = 0;
303676f2a66SJacob Faibussowitsch       if (++j <= right && ++gallopright >= MIN_GALLOP_GLOBAL) {
304676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
305676f2a66SJacob Faibussowitsch         do {
306676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
307676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
3084d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, j, right, atarr+(i*asize), &l2);CHKERRQ(ierr);
309676f2a66SJacob Faibussowitsch           diff2 = l2-j;
310676f2a66SJacob Faibussowitsch           Petsc_memmove2(arr+(k*asize), arr+(j*asize), diff2*asize, barr+(k*bsize), barr+(j*bsize), diff2*bsize);
311676f2a66SJacob Faibussowitsch           k += diff2;
312676f2a66SJacob Faibussowitsch           j = l2;
313676f2a66SJacob Faibussowitsch           /* search temp for right[j], can copy up to that of temp into arr immediately */
3144d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(atarr, asize, cmp, ctx, i, llen-1, arr+(j*asize), &l1);CHKERRQ(ierr);
315676f2a66SJacob Faibussowitsch           diff1 = l1-i;
316676f2a66SJacob Faibussowitsch           Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), diff1*asize, barr+(k*bsize), btarr+(i*bsize), diff1*bsize);
317676f2a66SJacob Faibussowitsch           k += diff1;
318676f2a66SJacob Faibussowitsch           i = l1;
319676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
320676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
321676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
322676f2a66SJacob Faibussowitsch       }
323676f2a66SJacob Faibussowitsch     }
324676f2a66SJacob Faibussowitsch   }
325676f2a66SJacob Faibussowitsch   if (i<llen) {Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), (llen-i)*asize, barr+(k*bsize), btarr+(i*bsize), (llen-i)*bsize);}
326676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
327676f2a66SJacob Faibussowitsch }
328676f2a66SJacob Faibussowitsch 
329676f2a66SJacob Faibussowitsch /* Mergesort where size of right half < size of left half, so mergesort is done right to left. Arr should be pointer to
330676f2a66SJacob Faibussowitsch  complete array, left is first index of left array, mid is first index of right array, right is last index of right
331676f2a66SJacob Faibussowitsch  array */
3324d3610e3SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortMergeHi_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt mid, PetscInt right)
333676f2a66SJacob Faibussowitsch {
334676f2a66SJacob Faibussowitsch   PetscInt       i = right-mid, j = mid-1, k = right, gallopleft = 0, gallopright = 0;
335676f2a66SJacob Faibussowitsch   const PetscInt rlen = right-mid+1;
336676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
337676f2a66SJacob Faibussowitsch 
338676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
339676f2a66SJacob Faibussowitsch   Petsc_memcpy(tarr, (arr)+mid*size, rlen*size);
340676f2a66SJacob Faibussowitsch   while ((i >= 0) && (j >= left)) {
3414d3610e3SJacob Faibussowitsch     if ((*cmp)((tarr)+i*size, (arr)+j*size, ctx) > 0) {
342676f2a66SJacob Faibussowitsch       Petsc_memcpy((arr)+k*size, (tarr)+i*size, size);
343676f2a66SJacob Faibussowitsch       --k;
344676f2a66SJacob Faibussowitsch       gallopleft = 0;
345676f2a66SJacob Faibussowitsch       if (--i >= 0 && ++gallopright >= MIN_GALLOP_GLOBAL) {
346676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
347676f2a66SJacob Faibussowitsch         do {
348676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
349676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
3504d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(tarr, size, cmp, ctx, 0, i, (arr)+j*size, &l1);CHKERRQ(ierr);
351676f2a66SJacob Faibussowitsch           diff1 = i-l1;
352676f2a66SJacob Faibussowitsch           Petsc_memcpy((arr)+(k-diff1+1)*size, (tarr)+(l1+1)*size, diff1*size);
353676f2a66SJacob Faibussowitsch           k -= diff1;
354676f2a66SJacob Faibussowitsch           i = l1;
355676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
3564d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, size, cmp, ctx, left, j, (tarr)+i*size, &l2);CHKERRQ(ierr);
357676f2a66SJacob Faibussowitsch           diff2 = j-l2;
358676f2a66SJacob Faibussowitsch           Petsc_memmove((arr)+(k-diff2+1)*size, (arr)+(l2+1)*size, diff2*size);
359676f2a66SJacob Faibussowitsch           k -= diff2;
360676f2a66SJacob Faibussowitsch           j = l2;
361676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
362676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
363676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
364676f2a66SJacob Faibussowitsch       }
365676f2a66SJacob Faibussowitsch     } else {
366676f2a66SJacob Faibussowitsch       Petsc_memmove((arr)+k*size, (arr)+j*size, size);
367676f2a66SJacob Faibussowitsch       --k;
368676f2a66SJacob Faibussowitsch       gallopright = 0;
369676f2a66SJacob Faibussowitsch       if (--j >= left && ++gallopleft >= MIN_GALLOP_GLOBAL) {
370676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
371676f2a66SJacob Faibussowitsch         do {
372676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
373676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
3744d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, size, cmp, ctx, left, j, (tarr)+i*size, &l2);CHKERRQ(ierr);
375676f2a66SJacob Faibussowitsch           diff2 = j-l2;
376676f2a66SJacob Faibussowitsch           Petsc_memmove((arr)+(k-diff2+1)*size, (arr)+(l2+1)*size, diff2*size);
377676f2a66SJacob Faibussowitsch           k -= diff2;
378676f2a66SJacob Faibussowitsch           j = l2;
379676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
3804d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(tarr, size, cmp, ctx, 0, i, (arr)+j*size, &l1);CHKERRQ(ierr);
381676f2a66SJacob Faibussowitsch           diff1 = i-l1;
382676f2a66SJacob Faibussowitsch           Petsc_memcpy((arr)+(k-diff1+1)*size, (tarr)+(l1+1)*size, diff1*size);
383676f2a66SJacob Faibussowitsch           k -= diff1;
384676f2a66SJacob Faibussowitsch           i = l1;
385676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
386676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
387676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
388676f2a66SJacob Faibussowitsch       }
389676f2a66SJacob Faibussowitsch     }
390676f2a66SJacob Faibussowitsch   }
391676f2a66SJacob Faibussowitsch   if (i >= 0) {Petsc_memcpy((arr)+left*size, tarr, (i+1)*size);}
392676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
393676f2a66SJacob Faibussowitsch }
394676f2a66SJacob Faibussowitsch 
3954d3610e3SJacob Faibussowitsch PETSC_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)
396676f2a66SJacob Faibussowitsch {
397676f2a66SJacob Faibussowitsch   PetscInt       i = right-mid, j = mid-1, k = right, gallopleft = 0, gallopright = 0;
398676f2a66SJacob Faibussowitsch   const PetscInt rlen = right-mid+1;
399676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
400676f2a66SJacob Faibussowitsch 
401676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
402676f2a66SJacob Faibussowitsch   Petsc_memcpy2(atarr, (arr)+mid*asize, rlen*asize, btarr, (barr)+mid*bsize, rlen*bsize);
403676f2a66SJacob Faibussowitsch   while ((i >= 0) && (j >= left)) {
4044d3610e3SJacob Faibussowitsch     if ((*cmp)((atarr)+i*asize, (arr)+j*asize, ctx) > 0) {
405676f2a66SJacob Faibussowitsch       Petsc_memcpy2((arr)+k*asize, (atarr)+i*asize, asize, (barr)+k*bsize, (btarr)+i*bsize, bsize);
406676f2a66SJacob Faibussowitsch       --k;
407676f2a66SJacob Faibussowitsch       gallopleft = 0;
408676f2a66SJacob Faibussowitsch       if (--i >= 0 && ++gallopright >= MIN_GALLOP_GLOBAL) {
409676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
410676f2a66SJacob Faibussowitsch         do {
411676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
412676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
4134d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(atarr, asize, cmp, ctx, 0, i, (arr)+j*asize, &l1);CHKERRQ(ierr);
414676f2a66SJacob Faibussowitsch           diff1 = i-l1;
415676f2a66SJacob 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);
416676f2a66SJacob Faibussowitsch           k -= diff1;
417676f2a66SJacob Faibussowitsch           i = l1;
418676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
4194d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, asize, cmp, ctx, left, j, (atarr)+i*asize, &l2);CHKERRQ(ierr);
420676f2a66SJacob Faibussowitsch           diff2 = j-l2;
421676f2a66SJacob 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);
422676f2a66SJacob Faibussowitsch           k -= diff2;
423676f2a66SJacob Faibussowitsch           j = l2;
424676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
425676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
426676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
427676f2a66SJacob Faibussowitsch       }
428676f2a66SJacob Faibussowitsch     } else {
429676f2a66SJacob Faibussowitsch       Petsc_memmove2((arr)+k*asize, (arr)+j*asize, asize, (barr)+k*bsize, (barr)+j*bsize, bsize);
430676f2a66SJacob Faibussowitsch       --k;
431676f2a66SJacob Faibussowitsch       gallopright = 0;
432676f2a66SJacob Faibussowitsch       if (--j >= left && ++gallopleft >= MIN_GALLOP_GLOBAL) {
433676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
434676f2a66SJacob Faibussowitsch         do {
435676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
436676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
4374d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, asize, cmp, ctx, left, j, (atarr)+i*asize, &l2);CHKERRQ(ierr);
438676f2a66SJacob Faibussowitsch           diff2 = j-l2;
439676f2a66SJacob 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);
440676f2a66SJacob Faibussowitsch           k -= diff2;
441676f2a66SJacob Faibussowitsch           j = l2;
442676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
4434d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(atarr, asize, cmp, ctx, 0, i, (arr)+j*asize, &l1);CHKERRQ(ierr);
444676f2a66SJacob Faibussowitsch           diff1 = i-l1;
445676f2a66SJacob 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);
446676f2a66SJacob Faibussowitsch           k -= diff1;
447676f2a66SJacob Faibussowitsch           i = l1;
448676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
449676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
450676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
451676f2a66SJacob Faibussowitsch       }
452676f2a66SJacob Faibussowitsch     }
453676f2a66SJacob Faibussowitsch   }
454676f2a66SJacob Faibussowitsch   if (i >= 0) {Petsc_memcpy2((arr)+left*asize, atarr, (i+1)*asize, (barr)+left*bsize, btarr, (i+1)*bsize);}
455676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
456676f2a66SJacob Faibussowitsch }
457676f2a66SJacob Faibussowitsch 
458676f2a66SJacob Faibussowitsch /* Left is inclusive lower bound of array slice, start is start location of unsorted section, right is inclusive upper
459676f2a66SJacob Faibussowitsch  bound of array slice. If unsure of where unsorted section starts or if entire length is unsorted pass start = left */
4604d3610e3SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscInsertionSort_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right)
461676f2a66SJacob Faibussowitsch {
462676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
463676f2a66SJacob Faibussowitsch 
464676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
465676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
466676f2a66SJacob Faibussowitsch     PetscInt j = i-1;
467676f2a66SJacob Faibussowitsch     Petsc_memcpy(tarr, arr+(i*size), size);
4684d3610e3SJacob Faibussowitsch     while ((j >= left) && ((*cmp)(tarr, (arr)+j*size, ctx) < 0)) {
469676f2a66SJacob Faibussowitsch       COPYSWAPPY(arr+(j+1)*size, arr+j*size, tarr+size, size);
470676f2a66SJacob Faibussowitsch       --j;
471676f2a66SJacob Faibussowitsch     }
472676f2a66SJacob Faibussowitsch     Petsc_memcpy((arr)+(j+1)*size, tarr, size);
473676f2a66SJacob Faibussowitsch   }
474676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
475676f2a66SJacob Faibussowitsch }
476676f2a66SJacob Faibussowitsch 
4774d3610e3SJacob Faibussowitsch PETSC_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)
478676f2a66SJacob Faibussowitsch {
479676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
480676f2a66SJacob Faibussowitsch 
481676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
482676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
483676f2a66SJacob Faibussowitsch     PetscInt j = i-1;
484676f2a66SJacob Faibussowitsch     Petsc_memcpy2(atarr, arr+(i*asize), asize, btarr, barr+(i*bsize), bsize);
4854d3610e3SJacob Faibussowitsch     while ((j >= left) && ((*cmp)(atarr, arr+(j*asize), ctx) < 0)) {
486676f2a66SJacob Faibussowitsch       COPYSWAPPY2(arr+(j+1)*asize, arr+j*asize, asize, barr+(j+1)*bsize, barr+j*bsize, bsize, atarr+asize);
487676f2a66SJacob Faibussowitsch       --j;
488676f2a66SJacob Faibussowitsch     }
489676f2a66SJacob Faibussowitsch     Petsc_memcpy2(arr+(j+1)*asize, atarr, asize, barr+(j+1)*bsize, btarr, bsize);
490676f2a66SJacob Faibussowitsch   }
491676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
492676f2a66SJacob Faibussowitsch }
493676f2a66SJacob Faibussowitsch 
494676f2a66SJacob Faibussowitsch /* See PetscInsertionSort_Private */
4954d3610e3SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscBinaryInsertionSort_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt left, PetscInt start, PetscInt right)
496676f2a66SJacob Faibussowitsch {
497676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
498676f2a66SJacob Faibussowitsch 
499676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
500676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
501676f2a66SJacob Faibussowitsch     PetscInt l = left, r = i;
502676f2a66SJacob Faibussowitsch     Petsc_memcpy(tarr, arr+(i*size), size);
503676f2a66SJacob Faibussowitsch     do {
504676f2a66SJacob Faibussowitsch       const PetscInt m = l + ((r - l) >> 1);
5054d3610e3SJacob Faibussowitsch       if ((*cmp)(tarr, arr+(m*size), ctx) < 0) {
506676f2a66SJacob Faibussowitsch         r = m;
507676f2a66SJacob Faibussowitsch       } else {
508676f2a66SJacob Faibussowitsch         l = m + 1;
509676f2a66SJacob Faibussowitsch       }
510676f2a66SJacob Faibussowitsch     } while (l < r);
511676f2a66SJacob Faibussowitsch     Petsc_memmove(arr+((l+1)*size), arr+(l*size), (i-l)*size);
512676f2a66SJacob Faibussowitsch     Petsc_memcpy(arr+(l*size), tarr, size);
513676f2a66SJacob Faibussowitsch   }
514676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
515676f2a66SJacob Faibussowitsch }
516676f2a66SJacob Faibussowitsch 
5174d3610e3SJacob Faibussowitsch PETSC_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)
518676f2a66SJacob Faibussowitsch {
519676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
520676f2a66SJacob Faibussowitsch 
521676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
522676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
523676f2a66SJacob Faibussowitsch     PetscInt l = left, r = i;
524676f2a66SJacob Faibussowitsch     Petsc_memcpy2(atarr, arr+(i*asize), asize, btarr, barr+(i*bsize), bsize);
525676f2a66SJacob Faibussowitsch     do {
526676f2a66SJacob Faibussowitsch       const PetscInt m = l + ((r - l) >> 1);
5274d3610e3SJacob Faibussowitsch       if ((*cmp)(atarr, arr+(m*asize), ctx) < 0) {
528676f2a66SJacob Faibussowitsch         r = m;
529676f2a66SJacob Faibussowitsch       } else {
530676f2a66SJacob Faibussowitsch         l = m + 1;
531676f2a66SJacob Faibussowitsch       }
532676f2a66SJacob Faibussowitsch     } while (l < r);
533676f2a66SJacob Faibussowitsch     Petsc_memmove2(arr+((l+1)*asize), arr+(l*asize), (i-l)*asize, barr+((l+1)*bsize), barr+(l*bsize), (i-l)*bsize);
534676f2a66SJacob Faibussowitsch     Petsc_memcpy2(arr+(l*asize), atarr, asize, barr+(l*bsize), btarr, bsize);
535676f2a66SJacob Faibussowitsch   }
536676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
537676f2a66SJacob Faibussowitsch }
538676f2a66SJacob Faibussowitsch 
539676f2a66SJacob Faibussowitsch typedef struct {
540676f2a66SJacob Faibussowitsch   PetscInt size;
541676f2a66SJacob Faibussowitsch   PetscInt start;
54289d292beSStefano Zampini #if defined(__CRAY_AARCH64) /* segfaults with Cray compilers for aarch64 on a64FX */
54389d292beSStefano Zampini } PetscTimSortStack;
54489d292beSStefano Zampini #else
545676f2a66SJacob Faibussowitsch } PetscTimSortStack PETSC_ATTRIBUTEALIGNED(2*sizeof(PetscInt));
54689d292beSStefano Zampini #endif
547676f2a66SJacob Faibussowitsch 
548676f2a66SJacob Faibussowitsch typedef struct {
549676f2a66SJacob Faibussowitsch   char   *ptr PETSC_ATTRIBUTEALIGNED(PETSC_MEMALIGN);
550676f2a66SJacob Faibussowitsch   size_t size;
551676f2a66SJacob Faibussowitsch   size_t maxsize;
552676f2a66SJacob Faibussowitsch } PetscTimSortBuffer;
553676f2a66SJacob Faibussowitsch 
554676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortResizeBuffer_Private(PetscTimSortBuffer *buff, size_t newSize)
555676f2a66SJacob Faibussowitsch {
556676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
557676f2a66SJacob Faibussowitsch   if (PetscLikely(newSize <= buff->size)) PetscFunctionReturn(0);
558676f2a66SJacob Faibussowitsch   {
559676f2a66SJacob Faibussowitsch     /* Can't be larger than n, there is merit to simply allocating buff to n to begin with */
560676f2a66SJacob Faibussowitsch     PetscErrorCode ierr;
561676f2a66SJacob Faibussowitsch     size_t         newMax = PetscMin(newSize*newSize, buff->maxsize);
562676f2a66SJacob Faibussowitsch     ierr = PetscFree(buff->ptr);CHKERRQ(ierr);
563676f2a66SJacob Faibussowitsch     ierr = PetscMalloc1(newMax, &buff->ptr);CHKERRQ(ierr);
564676f2a66SJacob Faibussowitsch     buff->size = newMax;
565676f2a66SJacob Faibussowitsch   }
566676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
567676f2a66SJacob Faibussowitsch }
568676f2a66SJacob Faibussowitsch 
5694d3610e3SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortForceCollapse_Private(char *arr, size_t size, CompFunc cmp, void *ctx, PetscTimSortBuffer *buff, PetscTimSortStack *stack, PetscInt stacksize)
570676f2a66SJacob Faibussowitsch {
571676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
572676f2a66SJacob Faibussowitsch   for (;stacksize; --stacksize) {
573676f2a66SJacob Faibussowitsch     /* A = stack[i-1], B = stack[i], if A[-1] <= B[0] means sorted */
5744d3610e3SJacob Faibussowitsch     if ((*cmp)(arr+(stack[stacksize].start-1)*size, arr+(stack[stacksize].start)*size, ctx) > 0) {
575676f2a66SJacob Faibussowitsch       PetscInt       l, m = stack[stacksize].start, r;
576676f2a66SJacob Faibussowitsch       PetscErrorCode ierr;
577676f2a66SJacob Faibussowitsch       /* Search A for B[0] insertion */
5784d3610e3SJacob Faibussowitsch       ierr = PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[stacksize-1].start, stack[stacksize].start-1, (arr)+(stack[stacksize].start)*size, &l);CHKERRQ(ierr);
579676f2a66SJacob Faibussowitsch       /* Search B for A[-1] insertion */
5804d3610e3SJacob Faibussowitsch       ierr = PetscGallopSearchRight_Private(arr, size, cmp, ctx, stack[stacksize].start, stack[stacksize].start+stack[stacksize].size-1, (arr)+(stack[stacksize].start-1)*size, &r);CHKERRQ(ierr);
581676f2a66SJacob Faibussowitsch       if (m-l <= r-m) {
582676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size);CHKERRQ(ierr);
5834d3610e3SJacob Faibussowitsch         ierr = PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r);CHKERRQ(ierr);
584676f2a66SJacob Faibussowitsch       } else {
585676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size);CHKERRQ(ierr);
5864d3610e3SJacob Faibussowitsch         ierr = PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r);CHKERRQ(ierr);
587676f2a66SJacob Faibussowitsch       }
588676f2a66SJacob Faibussowitsch     }
589676f2a66SJacob Faibussowitsch     /* Update A with merge */
590676f2a66SJacob Faibussowitsch     stack[stacksize-1].size += stack[stacksize].size;
591676f2a66SJacob Faibussowitsch   }
592676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
593676f2a66SJacob Faibussowitsch }
594676f2a66SJacob Faibussowitsch 
5954d3610e3SJacob Faibussowitsch PETSC_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)
596676f2a66SJacob Faibussowitsch {
597676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
598676f2a66SJacob Faibussowitsch   for (;stacksize; --stacksize) {
599676f2a66SJacob Faibussowitsch     /* A = stack[i-1], B = stack[i], if A[-1] <= B[0] means sorted */
6004d3610e3SJacob Faibussowitsch     if ((*cmp)(arr+(stack[stacksize].start-1)*asize, arr+(stack[stacksize].start)*asize, ctx) > 0) {
601676f2a66SJacob Faibussowitsch       PetscInt       l, m = stack[stacksize].start, r;
602676f2a66SJacob Faibussowitsch       PetscErrorCode ierr;
603676f2a66SJacob Faibussowitsch       /* Search A for B[0] insertion */
6044d3610e3SJacob Faibussowitsch       ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[stacksize-1].start, stack[stacksize].start-1, (arr)+(stack[stacksize].start)*asize, &l);CHKERRQ(ierr);
605676f2a66SJacob Faibussowitsch       /* Search B for A[-1] insertion */
6064d3610e3SJacob Faibussowitsch       ierr = PetscGallopSearchRight_Private(arr, asize, cmp, ctx, stack[stacksize].start, stack[stacksize].start+stack[stacksize].size-1, (arr)+(stack[stacksize].start-1)*asize, &r);CHKERRQ(ierr);
607676f2a66SJacob Faibussowitsch       if (m-l <= r-m) {
608676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize);CHKERRQ(ierr);
609676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize);CHKERRQ(ierr);
6104d3610e3SJacob Faibussowitsch         ierr = PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r);CHKERRQ(ierr);
611676f2a66SJacob Faibussowitsch       } else {
612676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize);CHKERRQ(ierr);
613676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize);CHKERRQ(ierr);
6144d3610e3SJacob Faibussowitsch         ierr = PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r);CHKERRQ(ierr);
615676f2a66SJacob Faibussowitsch       }
616676f2a66SJacob Faibussowitsch     }
617676f2a66SJacob Faibussowitsch     /* Update A with merge */
618676f2a66SJacob Faibussowitsch     stack[stacksize-1].size += stack[stacksize].size;
619676f2a66SJacob Faibussowitsch   }
620676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
621676f2a66SJacob Faibussowitsch }
622676f2a66SJacob Faibussowitsch 
6234d3610e3SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortMergeCollapse_Private(char *arr, size_t size, CompFunc cmp, void *ctx, PetscTimSortBuffer *buff, PetscTimSortStack *stack, PetscInt *stacksize)
624676f2a66SJacob Faibussowitsch {
625676f2a66SJacob Faibussowitsch   PetscInt       i = *stacksize;
626676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
627676f2a66SJacob Faibussowitsch 
628676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
629676f2a66SJacob Faibussowitsch   while (i) {
630676f2a66SJacob Faibussowitsch     PetscInt l, m, r, itemp = i;
631676f2a66SJacob Faibussowitsch 
632676f2a66SJacob Faibussowitsch     if (i == 1) {
633676f2a66SJacob Faibussowitsch       /* A = stack[i-1], B = stack[i] */
634676f2a66SJacob Faibussowitsch       if (stack[i-1].size < stack[i].size) {
635676f2a66SJacob Faibussowitsch         /* if A[-1] <= B[0] then sorted */
6364d3610e3SJacob Faibussowitsch         if ((*cmp)(arr+(stack[i].start-1)*size, arr+(stack[i].start)*size, ctx) > 0) {
637676f2a66SJacob Faibussowitsch           m = stack[i].start;
638676f2a66SJacob Faibussowitsch           /* Search A for B[0] insertion */
6394d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i-1].start, stack[i].start-1, (arr)+stack[i].start*size, &l);CHKERRQ(ierr);
640676f2a66SJacob Faibussowitsch           /* Search B for A[-1] insertion */
6414d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, size, cmp, ctx, stack[i].start, stack[i].start+stack[i].size-1, arr+(stack[i].start-1)*size, &r);CHKERRQ(ierr);
642676f2a66SJacob Faibussowitsch           if (m-l <= r-m) {
643676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size);CHKERRQ(ierr);
6444d3610e3SJacob Faibussowitsch             ierr = PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r);CHKERRQ(ierr);
645676f2a66SJacob Faibussowitsch           } else {
646676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size);CHKERRQ(ierr);
6474d3610e3SJacob Faibussowitsch             ierr = PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r);CHKERRQ(ierr);
648676f2a66SJacob Faibussowitsch           }
649676f2a66SJacob Faibussowitsch         }
650676f2a66SJacob Faibussowitsch         /* Update A with merge */
651676f2a66SJacob Faibussowitsch         stack[i-1].size += stack[i].size;
652676f2a66SJacob Faibussowitsch         --i;
653676f2a66SJacob Faibussowitsch       }
654676f2a66SJacob Faibussowitsch     } else {
655676f2a66SJacob Faibussowitsch       /* i > 2, i.e. C exists
656676f2a66SJacob Faibussowitsch        A = stack[i-2], B = stack[i-1], C = stack[i]; */
657676f2a66SJacob Faibussowitsch       if (stack[i-2].size <= stack[i-1].size+stack[i].size) {
658676f2a66SJacob Faibussowitsch         if (stack[i-2].size < stack[i].size) {
659676f2a66SJacob Faibussowitsch           /* merge B into A, but if A[-1] <= B[0] then already sorted */
6604d3610e3SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i-1].start-1)*size, arr+(stack[i-1].start)*size, ctx) > 0) {
661676f2a66SJacob Faibussowitsch             m = stack[i-1].start;
662676f2a66SJacob Faibussowitsch             /* Search A for B[0] insertion */
6634d3610e3SJacob Faibussowitsch             ierr = PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i-2].start, stack[i-1].start-1, (arr)+(stack[i-1].start)*size, &l);CHKERRQ(ierr);
664676f2a66SJacob Faibussowitsch             /* Search B for A[-1] insertion */
6654d3610e3SJacob Faibussowitsch             ierr = 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);CHKERRQ(ierr);
666676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
667676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size);CHKERRQ(ierr);
6684d3610e3SJacob Faibussowitsch               ierr = PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r);CHKERRQ(ierr);
669676f2a66SJacob Faibussowitsch             } else {
670676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size);CHKERRQ(ierr);
6714d3610e3SJacob Faibussowitsch               ierr = PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r);CHKERRQ(ierr);
672676f2a66SJacob Faibussowitsch             }
673676f2a66SJacob Faibussowitsch           }
674676f2a66SJacob Faibussowitsch           /* Update A with merge */
675676f2a66SJacob Faibussowitsch           stack[i-2].size += stack[i-1].size;
676676f2a66SJacob Faibussowitsch           /* Push C up the stack */
677676f2a66SJacob Faibussowitsch           stack[i-1].start = stack[i].start;
678676f2a66SJacob Faibussowitsch           stack[i-1].size = stack[i].size;
679676f2a66SJacob Faibussowitsch         } else {
680676f2a66SJacob Faibussowitsch           /* merge C into B */
681676f2a66SJacob Faibussowitsch           mergeBC:
682676f2a66SJacob Faibussowitsch           /* If B[-1] <= C[0] then... you know the drill */
6834d3610e3SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i].start-1)*size, arr+(stack[i].start)*size, ctx) > 0) {
684676f2a66SJacob Faibussowitsch             m = stack[i].start;
685676f2a66SJacob Faibussowitsch             /* Search B for C[0] insertion */
6864d3610e3SJacob Faibussowitsch             ierr = PetscGallopSearchLeft_Private(arr, size, cmp, ctx, stack[i-1].start, stack[i].start-1, arr+stack[i].start*size, &l);CHKERRQ(ierr);
687676f2a66SJacob Faibussowitsch             /* Search C for B[-1] insertion */
6884d3610e3SJacob Faibussowitsch             ierr = PetscGallopSearchRight_Private(arr, size, cmp, ctx, stack[i].start, stack[i].start+stack[i].size-1, (arr)+(stack[i].start-1)*size, &r);CHKERRQ(ierr);
689676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
690676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size);CHKERRQ(ierr);
6914d3610e3SJacob Faibussowitsch               ierr = PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, ctx, l, m, r);CHKERRQ(ierr);
692676f2a66SJacob Faibussowitsch             } else {
693676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size);CHKERRQ(ierr);
6944d3610e3SJacob Faibussowitsch               ierr = PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, ctx, l, m, r);CHKERRQ(ierr);
695676f2a66SJacob Faibussowitsch             }
696676f2a66SJacob Faibussowitsch           }
697676f2a66SJacob Faibussowitsch           /* Update B with merge */
698676f2a66SJacob Faibussowitsch           stack[i-1].size += stack[i].size;
699676f2a66SJacob Faibussowitsch         }
700676f2a66SJacob Faibussowitsch         --i;
701676f2a66SJacob Faibussowitsch       } else if (stack[i-1].size <= stack[i].size) {
702676f2a66SJacob Faibussowitsch         /* merge C into B */
703676f2a66SJacob Faibussowitsch         goto mergeBC;
704676f2a66SJacob Faibussowitsch       }
705676f2a66SJacob Faibussowitsch     }
706676f2a66SJacob Faibussowitsch     if (itemp == i) break;
707676f2a66SJacob Faibussowitsch   }
708676f2a66SJacob Faibussowitsch   *stacksize = i;
709676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
710676f2a66SJacob Faibussowitsch }
711676f2a66SJacob Faibussowitsch 
7124d3610e3SJacob Faibussowitsch PETSC_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)
713676f2a66SJacob Faibussowitsch {
714676f2a66SJacob Faibussowitsch   PetscInt       i = *stacksize;
715676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
716676f2a66SJacob Faibussowitsch 
717676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
718676f2a66SJacob Faibussowitsch   while (i) {
719676f2a66SJacob Faibussowitsch     PetscInt l, m, r, itemp = i;
720676f2a66SJacob Faibussowitsch 
721676f2a66SJacob Faibussowitsch     if (i == 1) {
722676f2a66SJacob Faibussowitsch       /* A = stack[i-1], B = stack[i] */
723676f2a66SJacob Faibussowitsch       if (stack[i-1].size < stack[i].size) {
724676f2a66SJacob Faibussowitsch         /* if A[-1] <= B[0] then sorted */
7254d3610e3SJacob Faibussowitsch         if ((*cmp)(arr+(stack[i].start-1)*asize, arr+(stack[i].start)*asize, ctx) > 0) {
726676f2a66SJacob Faibussowitsch           m = stack[i].start;
727676f2a66SJacob Faibussowitsch           /* Search A for B[0] insertion */
7284d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i-1].start, stack[i].start-1, (arr)+stack[i].start*asize, &l);CHKERRQ(ierr);
729676f2a66SJacob Faibussowitsch           /* Search B for A[-1] insertion */
7304d3610e3SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, asize, cmp, ctx, stack[i].start, stack[i].start+stack[i].size-1, arr+(stack[i].start-1)*asize, &r);CHKERRQ(ierr);
731676f2a66SJacob Faibussowitsch           if (m-l <= r-m) {
732676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize);CHKERRQ(ierr);
733676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize);CHKERRQ(ierr);
7344d3610e3SJacob Faibussowitsch             ierr = PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r);CHKERRQ(ierr);
735676f2a66SJacob Faibussowitsch           } else {
736676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize);CHKERRQ(ierr);
737676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize);CHKERRQ(ierr);
7384d3610e3SJacob Faibussowitsch             ierr = PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r);CHKERRQ(ierr);
739676f2a66SJacob Faibussowitsch           }
740676f2a66SJacob Faibussowitsch         }
741676f2a66SJacob Faibussowitsch         /* Update A with merge */
742676f2a66SJacob Faibussowitsch         stack[i-1].size += stack[i].size;
743676f2a66SJacob Faibussowitsch         --i;
744676f2a66SJacob Faibussowitsch       }
745676f2a66SJacob Faibussowitsch     } else {
746676f2a66SJacob Faibussowitsch       /* i > 2, i.e. C exists
747676f2a66SJacob Faibussowitsch        A = stack[i-2], B = stack[i-1], C = stack[i]; */
748676f2a66SJacob Faibussowitsch       if (stack[i-2].size <= stack[i-1].size+stack[i].size) {
749676f2a66SJacob Faibussowitsch         if (stack[i-2].size < stack[i].size) {
750676f2a66SJacob Faibussowitsch           /* merge B into A, but if A[-1] <= B[0] then already sorted */
7514d3610e3SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i-1].start-1)*asize, arr+(stack[i-1].start)*asize, ctx) > 0) {
752676f2a66SJacob Faibussowitsch             m = stack[i-1].start;
753676f2a66SJacob Faibussowitsch             /* Search A for B[0] insertion */
7544d3610e3SJacob Faibussowitsch             ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i-2].start, stack[i-1].start-1, (arr)+(stack[i-1].start)*asize, &l);CHKERRQ(ierr);
755676f2a66SJacob Faibussowitsch             /* Search B for A[-1] insertion */
7564d3610e3SJacob Faibussowitsch             ierr = 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);CHKERRQ(ierr);
757676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
758676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize);CHKERRQ(ierr);
759676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize);CHKERRQ(ierr);
7604d3610e3SJacob Faibussowitsch               ierr = PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r);CHKERRQ(ierr);
761676f2a66SJacob Faibussowitsch             } else {
762676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize);CHKERRQ(ierr);
763676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize);CHKERRQ(ierr);
7644d3610e3SJacob Faibussowitsch               ierr = PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r);CHKERRQ(ierr);
765676f2a66SJacob Faibussowitsch             }
766676f2a66SJacob Faibussowitsch           }
767676f2a66SJacob Faibussowitsch           /* Update A with merge */
768676f2a66SJacob Faibussowitsch           stack[i-2].size += stack[i-1].size;
769676f2a66SJacob Faibussowitsch           /* Push C up the stack */
770676f2a66SJacob Faibussowitsch           stack[i-1].start = stack[i].start;
771676f2a66SJacob Faibussowitsch           stack[i-1].size = stack[i].size;
772676f2a66SJacob Faibussowitsch         } else {
773676f2a66SJacob Faibussowitsch           /* merge C into B */
774676f2a66SJacob Faibussowitsch           mergeBC:
775676f2a66SJacob Faibussowitsch           /* If B[-1] <= C[0] then... you know the drill */
7764d3610e3SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i].start-1)*asize, arr+(stack[i].start)*asize, ctx) > 0) {
777676f2a66SJacob Faibussowitsch             m = stack[i].start;
778676f2a66SJacob Faibussowitsch             /* Search B for C[0] insertion */
7794d3610e3SJacob Faibussowitsch             ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, ctx, stack[i-1].start, stack[i].start-1, arr+stack[i].start*asize, &l);CHKERRQ(ierr);
780676f2a66SJacob Faibussowitsch             /* Search C for B[-1] insertion */
7814d3610e3SJacob Faibussowitsch             ierr = PetscGallopSearchRight_Private(arr, asize, cmp, ctx, stack[i].start, stack[i].start+stack[i].size-1, (arr)+(stack[i].start-1)*asize, &r);CHKERRQ(ierr);
782676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
783676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize);CHKERRQ(ierr);
784676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize);CHKERRQ(ierr);
7854d3610e3SJacob Faibussowitsch               ierr = PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r);CHKERRQ(ierr);
786676f2a66SJacob Faibussowitsch             } else {
787676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize);CHKERRQ(ierr);
788676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize);CHKERRQ(ierr);
7894d3610e3SJacob Faibussowitsch               ierr = PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, ctx, l, m, r);CHKERRQ(ierr);
790676f2a66SJacob Faibussowitsch             }
791676f2a66SJacob Faibussowitsch           }
792676f2a66SJacob Faibussowitsch           /* Update B with merge */
793676f2a66SJacob Faibussowitsch           stack[i-1].size += stack[i].size;
794676f2a66SJacob Faibussowitsch         }
795676f2a66SJacob Faibussowitsch         --i;
796676f2a66SJacob Faibussowitsch       } else if (stack[i-1].size <= stack[i].size) {
797676f2a66SJacob Faibussowitsch         /* merge C into B */
798676f2a66SJacob Faibussowitsch         goto mergeBC;
799676f2a66SJacob Faibussowitsch       }
800676f2a66SJacob Faibussowitsch     }
801676f2a66SJacob Faibussowitsch     if (itemp == i) break;
802676f2a66SJacob Faibussowitsch   }
803676f2a66SJacob Faibussowitsch   *stacksize = i;
804676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
805676f2a66SJacob Faibussowitsch }
806676f2a66SJacob Faibussowitsch 
807676f2a66SJacob Faibussowitsch /* March sequentially through the array building up a "run" of weakly increasing or strictly decreasing contiguous
808676f2a66SJacob Faibussowitsch  elements. Decreasing runs are reversed by swapping. If the run is less than minrun, artificially extend it via either
809676f2a66SJacob Faibussowitsch  binary insertion sort or regulat insertion sort */
8104d3610e3SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortBuildRun_Private(char *arr, char *tarr, size_t size, CompFunc cmp, void *ctx, PetscInt n, PetscInt minrun, PetscInt runstart, PetscInt *runend)
811676f2a66SJacob Faibussowitsch {
812676f2a66SJacob Faibussowitsch   const PetscInt re = PetscMin(runstart+minrun, n-1);
813676f2a66SJacob Faibussowitsch   PetscInt       ri = runstart;
814676f2a66SJacob Faibussowitsch 
815676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
816676f2a66SJacob Faibussowitsch   if (PetscUnlikely(runstart == n-1)) {*runend = runstart; PetscFunctionReturn(0);}
817676f2a66SJacob Faibussowitsch   /* guess whether run is ascending or descending and tally up the longest consecutive run. essentially a coinflip for random data */
8184d3610e3SJacob Faibussowitsch   if ((*cmp)((arr)+(ri+1)*size, (arr)+ri*size, ctx) < 0) {
819676f2a66SJacob Faibussowitsch     ++ri;
820676f2a66SJacob Faibussowitsch     while (ri < n-1) {
8214d3610e3SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*size, (arr)+ri*size, ctx) >= 0) break;
822676f2a66SJacob Faibussowitsch       ++ri;
823676f2a66SJacob Faibussowitsch     }
824676f2a66SJacob Faibussowitsch     {
825676f2a66SJacob Faibussowitsch       PetscInt lo = runstart, hi = ri;
826676f2a66SJacob Faibussowitsch       for (; lo < hi; ++lo, --hi) {
827676f2a66SJacob Faibussowitsch         COPYSWAPPY(arr+lo*size, arr+hi*size, tarr, size);
828676f2a66SJacob Faibussowitsch       }
829676f2a66SJacob Faibussowitsch     }
830676f2a66SJacob Faibussowitsch   } else {
831676f2a66SJacob Faibussowitsch     ++ri;
832676f2a66SJacob Faibussowitsch     while (ri < n-1) {
8334d3610e3SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*size, (arr)+ri*size, ctx) < 0) break;
834676f2a66SJacob Faibussowitsch       ++ri;
835676f2a66SJacob Faibussowitsch     }
836676f2a66SJacob Faibussowitsch   }
837676f2a66SJacob Faibussowitsch #if defined(PETSC_USE_DEBUG)
838676f2a66SJacob Faibussowitsch   {
839676f2a66SJacob Faibussowitsch     PetscErrorCode ierr;
840*3ca90d2dSJacob Faibussowitsch     ierr = PetscInfo1(NULL, "natural run length = %" PetscInt_FMT "\n", ri-runstart+1);CHKERRQ(ierr);
841676f2a66SJacob Faibussowitsch   }
842676f2a66SJacob Faibussowitsch #endif
843676f2a66SJacob Faibussowitsch   if (ri < re) {
844676f2a66SJacob Faibussowitsch     /* the attempt failed, this section likely contains random data. If ri got close to minrun (within 50%) then we try
845676f2a66SJacob Faibussowitsch      binary search */
846676f2a66SJacob Faibussowitsch     if (ri-runstart <= minrun >> 1) {
847676f2a66SJacob Faibussowitsch       ++MIN_GALLOP_GLOBAL; /* didn't get close hedge our bets against random data */
8484d3610e3SJacob Faibussowitsch       PetscInsertionSort_Private(arr, tarr, size, cmp, ctx, runstart, ri, re);
849676f2a66SJacob Faibussowitsch     } else {
8504d3610e3SJacob Faibussowitsch       PetscBinaryInsertionSort_Private(arr, tarr, size, cmp, ctx, runstart, ri, re);
851676f2a66SJacob Faibussowitsch     }
852676f2a66SJacob Faibussowitsch     *runend = re;
853676f2a66SJacob Faibussowitsch   } else *runend = ri;
854676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
855676f2a66SJacob Faibussowitsch }
856676f2a66SJacob Faibussowitsch 
8574d3610e3SJacob Faibussowitsch PETSC_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)
858676f2a66SJacob Faibussowitsch {
859676f2a66SJacob Faibussowitsch   const PetscInt re = PetscMin(runstart+minrun, n-1);
860676f2a66SJacob Faibussowitsch   PetscInt       ri = runstart;
861676f2a66SJacob Faibussowitsch 
862676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
863676f2a66SJacob Faibussowitsch   if (PetscUnlikely(runstart == n-1)) {*runend = runstart; PetscFunctionReturn(0);}
864676f2a66SJacob Faibussowitsch   /* guess whether run is ascending or descending and tally up the longest consecutive run. essentially a coinflip for random data */
8654d3610e3SJacob Faibussowitsch   if ((*cmp)((arr)+(ri+1)*asize, arr+(ri*asize), ctx) < 0) {
866676f2a66SJacob Faibussowitsch     ++ri;
867676f2a66SJacob Faibussowitsch     while (ri < n-1) {
8684d3610e3SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*asize, (arr)+ri*asize, ctx) >= 0) break;
869676f2a66SJacob Faibussowitsch       ++ri;
870676f2a66SJacob Faibussowitsch     }
871676f2a66SJacob Faibussowitsch     {
872676f2a66SJacob Faibussowitsch       PetscInt lo = runstart, hi = ri;
873676f2a66SJacob Faibussowitsch       for (; lo < hi; ++lo, --hi) {
874676f2a66SJacob Faibussowitsch         COPYSWAPPY2(arr+lo*asize, arr+hi*asize, asize, barr+lo*bsize, barr+hi*bsize, bsize, atarr);
875676f2a66SJacob Faibussowitsch       }
876676f2a66SJacob Faibussowitsch     }
877676f2a66SJacob Faibussowitsch   } else {
878676f2a66SJacob Faibussowitsch     ++ri;
879676f2a66SJacob Faibussowitsch     while (ri < n-1) {
8804d3610e3SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*asize, (arr)+ri*asize, ctx) < 0) break;
881676f2a66SJacob Faibussowitsch       ++ri;
882676f2a66SJacob Faibussowitsch     }
883676f2a66SJacob Faibussowitsch   }
884676f2a66SJacob Faibussowitsch #if defined(PETSC_USE_DEBUG)
885676f2a66SJacob Faibussowitsch   {
886676f2a66SJacob Faibussowitsch     PetscErrorCode ierr;
887*3ca90d2dSJacob Faibussowitsch     ierr = PetscInfo1(NULL, "natural run length = %" PetscInt_FMT "\n", ri-runstart+1);CHKERRQ(ierr);
888676f2a66SJacob Faibussowitsch   }
889676f2a66SJacob Faibussowitsch #endif
890676f2a66SJacob Faibussowitsch   if (ri < re) {
891676f2a66SJacob Faibussowitsch     /* the attempt failed, this section likely contains random data. If ri got close to minrun (within 50%) then we try
892676f2a66SJacob Faibussowitsch      binary search */
893676f2a66SJacob Faibussowitsch     if (ri-runstart <= minrun >> 1) {
894676f2a66SJacob Faibussowitsch       ++MIN_GALLOP_GLOBAL; /* didn't get close hedge our bets against random data */
8954d3610e3SJacob Faibussowitsch       PetscInsertionSortWithArray_Private(arr, atarr, asize, barr, btarr, bsize, cmp, ctx, runstart, ri, re);
896676f2a66SJacob Faibussowitsch     } else {
8974d3610e3SJacob Faibussowitsch       PetscBinaryInsertionSortWithArray_Private(arr, atarr, asize, barr, btarr, bsize, cmp, ctx, runstart, ri, re);
898676f2a66SJacob Faibussowitsch     }
899676f2a66SJacob Faibussowitsch     *runend = re;
900676f2a66SJacob Faibussowitsch   } else *runend = ri;
901676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
902676f2a66SJacob Faibussowitsch }
903676f2a66SJacob Faibussowitsch 
9044d3610e3SJacob Faibussowitsch /*@C
905676f2a66SJacob Faibussowitsch   PetscTimSort - Sorts an array in place in increasing order using Tim Peters adaptive sorting algorithm.
906676f2a66SJacob Faibussowitsch 
907676f2a66SJacob Faibussowitsch   Not Collective
908676f2a66SJacob Faibussowitsch 
909676f2a66SJacob Faibussowitsch   Input Parameters:
910676f2a66SJacob Faibussowitsch + n    - number of values
911676f2a66SJacob Faibussowitsch . arr  - array to be sorted
912676f2a66SJacob Faibussowitsch . size - size in bytes of the datatype held in arr
9134d3610e3SJacob Faibussowitsch . cmp  - function pointer to comparison function
9144d3610e3SJacob Faibussowitsch - ctx  - optional context to be passed to comparison function, NULL if not needed
915676f2a66SJacob Faibussowitsch 
916676f2a66SJacob Faibussowitsch   Output Parameters:
917676f2a66SJacob Faibussowitsch . arr  - sorted array
918676f2a66SJacob Faibussowitsch 
9198da36234SJacob Faibussowitsch   Notes:
9208da36234SJacob Faibussowitsch   Timsort makes the assumption that input data is already likely partially ordered, or that it contains contiguous
9218da36234SJacob Faibussowitsch  sections (termed 'runs') where the data is locally ordered (but not necessarily globally ordered). It therefore aims to
9228da36234SJacob Faibussowitsch  select slices of the array in such a way that resulting mergesorts operate on near perfectly length-balanced arrays. To
9238da36234SJacob Faibussowitsch  do so it repeatedly triggers attempts throughout to merge adjacent runs.
9248da36234SJacob Faibussowitsch 
9258da36234SJacob Faibussowitsch   Should one run continuously "win" a comparison the algorithm begins the "gallop" phase. It will aggressively search
9268da36234SJacob Faibussowitsch   the "winner" for the location of the "losers" next entry (and vice versa) to copy all preceding elements into place in
9278da36234SJacob Faibussowitsch   bulk. However if the data is truly unordered (as is the case with random data) the immense gains possible from these
9288da36234SJacob Faibussowitsch   searches are expected __not__ to repay their costs. While adjacent arrays are almost all nearly the same size, they
9298da36234SJacob Faibussowitsch   likely all contain similar data.
9308da36234SJacob Faibussowitsch 
931676f2a66SJacob Faibussowitsch   Sample usage:
932676f2a66SJacob Faibussowitsch   The comparison function must follow the qsort() comparison function paradigm, returning the sign of the difference
933676f2a66SJacob Faibussowitsch   between its arguments. If left < right : return -1, if left == right : return 0, if left > right : return 1. The user
934676f2a66SJacob Faibussowitsch   may also
935676f2a66SJacob Faibussowitsch  change or reverse the order of the sort by flipping the above. Note that stability of the sort is only guaranteed if
936676f2a66SJacob Faibussowitsch  the comparison function forms a valid trigraph. For example when sorting an array of type "my_type" in increasing
9378da36234SJacob Faibussowitsch   order
938676f2a66SJacob Faibussowitsch 
939676f2a66SJacob Faibussowitsch .vb
9404d3610e3SJacob Faibussowitsch   int my_increasing_comparison_function(const void *left, const void *right, void *ctx) {
941676f2a66SJacob Faibussowitsch     my_type l = *(my_type *) left, r = *(my_type *) right;
942676f2a66SJacob Faibussowitsch     return (l < r) ? -1 : (l > r);
943676f2a66SJacob Faibussowitsch   }
944676f2a66SJacob Faibussowitsch .ve
9454d3610e3SJacob Faibussowitsch   Note the context is unused here but you may use it to pass and subsequently access whatever information required
9464d3610e3SJacob Faibussowitsch   inside the comparison function. The context pointer will unaltered except for any changes made inside the comparison function.
947676f2a66SJacob Faibussowitsch   Then pass the function
948676f2a66SJacob Faibussowitsch .vb
9494d3610e3SJacob Faibussowitsch   PetscTimSort(n, arr, sizeof(arr[0]), my_increasing_comparison_function, ctx)
950676f2a66SJacob Faibussowitsch .ve
951676f2a66SJacob Faibussowitsch 
9524d3610e3SJacob Faibussowitsch   Fortran Notes:
9534d3610e3SJacob Faibussowitsch   To use this from fortran you must write a comparison subroutine with 4 arguments which accepts left, right, ctx and
9548da36234SJacob Faibussowitsch   returns result. For example
9554d3610e3SJacob Faibussowitsch .vb
9564d3610e3SJacob Faibussowitsch  subroutine CompareIntegers(left,right,ctx,result)
9574d3610e3SJacob Faibussowitsch    implicit none
9584d3610e3SJacob Faibussowitsch 
9594d3610e3SJacob Faibussowitsch    PetscInt,intent(in) :: left, right
9604d3610e3SJacob Faibussowitsch    type(UserCtx)       :: ctx
9614d3610e3SJacob Faibussowitsch    integer,intent(out) :: result
9624d3610e3SJacob Faibussowitsch 
9634d3610e3SJacob Faibussowitsch    if (left < right) then
9644d3610e3SJacob Faibussowitsch      result = -1
9654d3610e3SJacob Faibussowitsch    else if (left == right) then
9664d3610e3SJacob Faibussowitsch      result = 0
9674d3610e3SJacob Faibussowitsch    else
9684d3610e3SJacob Faibussowitsch      result = 1
9694d3610e3SJacob Faibussowitsch    end if
9704d3610e3SJacob Faibussowitsch    return
9714d3610e3SJacob Faibussowitsch  end subroutine CompareIntegers
9724d3610e3SJacob Faibussowitsch .ve
9734d3610e3SJacob Faibussowitsch 
974676f2a66SJacob Faibussowitsch   References:
975676f2a66SJacob Faibussowitsch   1. - Tim Peters. https://bugs.python.org/file4451/timsort.txt
976676f2a66SJacob Faibussowitsch 
977676f2a66SJacob Faibussowitsch   Level: developer
978676f2a66SJacob Faibussowitsch 
979676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscIntSortSemiOrdered(), PetscRealSortSemiOrdered(), PetscMPIIntSortSemiOrdered()
9804d3610e3SJacob Faibussowitsch @*/
9814d3610e3SJacob Faibussowitsch PetscErrorCode PetscTimSort(PetscInt n, void *arr, size_t size, int (*cmp)(const void *, const void *, void *), void *ctx)
982676f2a66SJacob Faibussowitsch {
983676f2a66SJacob Faibussowitsch   PetscInt           stacksize = 0, minrun, runstart = 0, runend = 0;
984676f2a66SJacob Faibussowitsch   PetscTimSortStack  runstack[128];
985676f2a66SJacob Faibussowitsch   PetscTimSortBuffer buff;
986676f2a66SJacob Faibussowitsch   PetscErrorCode     ierr;
987676f2a66SJacob Faibussowitsch   /* stacksize  = log_phi(n) = log_2(n)/log_2(phi), so 128 is enough for ~5.614e26 elements.
988676f2a66SJacob Faibussowitsch    It is so unlikely that this limit is reached that this is __never__ checked for */
989676f2a66SJacob Faibussowitsch 
990676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
991676f2a66SJacob Faibussowitsch   /* Compute minrun. Minrun should be (32, 65) such that N/minrun
992676f2a66SJacob Faibussowitsch    is a power of 2 or one plus a power of 2 */
993676f2a66SJacob Faibussowitsch   {
994676f2a66SJacob Faibussowitsch     PetscInt t = n, r = 0;
995676f2a66SJacob Faibussowitsch     /* r becomes 1 if the least significant bits contain at least one off bit */
996676f2a66SJacob Faibussowitsch     while (t >= 64) {
997676f2a66SJacob Faibussowitsch       r |= t & 1;
998676f2a66SJacob Faibussowitsch       t >>= 1;
999676f2a66SJacob Faibussowitsch     }
1000676f2a66SJacob Faibussowitsch     minrun = t + r;
1001676f2a66SJacob Faibussowitsch   }
10024d3610e3SJacob Faibussowitsch   if (PetscDefined(USE_DEBUG)) {
1003*3ca90d2dSJacob Faibussowitsch     ierr = PetscInfo1(NULL, "minrun = %" PetscInt_FMT "\n", minrun);CHKERRQ(ierr);
10044d3610e3SJacob Faibussowitsch     if (n < 64) {
1005*3ca90d2dSJacob Faibussowitsch       ierr = PetscInfo1(NULL, "n %" PetscInt_FMT " < 64, consider using PetscSortInt() instead\n", n);CHKERRQ(ierr);
1006*3ca90d2dSJacob Faibussowitsch     } else if ((minrun < 32) || (minrun > 65)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Calculated minrun %" PetscInt_FMT " not in range (32,65)",minrun);
10074d3610e3SJacob Faibussowitsch   }
1008676f2a66SJacob Faibussowitsch   ierr = PetscMalloc1((size_t) minrun*size, &buff.ptr);CHKERRQ(ierr);
1009676f2a66SJacob Faibussowitsch   buff.size = (size_t) minrun*size;
1010676f2a66SJacob Faibussowitsch   buff.maxsize = (size_t) n*size;
1011676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
1012676f2a66SJacob Faibussowitsch   while (runstart < n) {
1013676f2a66SJacob Faibussowitsch     /* Check if additional entries are at least partially ordered and build natural run */
10144d3610e3SJacob Faibussowitsch     ierr = PetscTimSortBuildRun_Private((char *)arr, buff.ptr, size, cmp, ctx, n, minrun, runstart, &runend);CHKERRQ(ierr);
1015676f2a66SJacob Faibussowitsch     runstack[stacksize].start = runstart;
1016676f2a66SJacob Faibussowitsch     runstack[stacksize].size = runend-runstart+1;
10174d3610e3SJacob Faibussowitsch     ierr = PetscTimSortMergeCollapse_Private((char *)arr, size, cmp, ctx, &buff, runstack, &stacksize);CHKERRQ(ierr);
1018676f2a66SJacob Faibussowitsch     ++stacksize;
1019676f2a66SJacob Faibussowitsch     runstart = runend+1;
1020676f2a66SJacob Faibussowitsch   }
1021676f2a66SJacob Faibussowitsch   /* Have been inside while, so discard last stacksize++ */
1022676f2a66SJacob Faibussowitsch   --stacksize;
10234d3610e3SJacob Faibussowitsch   ierr = PetscTimSortForceCollapse_Private((char *)arr, size, cmp, ctx, &buff, runstack, stacksize);CHKERRQ(ierr);
1024676f2a66SJacob Faibussowitsch   ierr = PetscFree(buff.ptr);CHKERRQ(ierr);
1025676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
1026676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1027676f2a66SJacob Faibussowitsch }
1028676f2a66SJacob Faibussowitsch 
10294d3610e3SJacob Faibussowitsch /*@C
1030676f2a66SJacob Faibussowitsch   PetscTimSortWithArray - Sorts an array in place in increasing order using Tim Peters adaptive sorting algorithm and
1031676f2a66SJacob Faibussowitsch   reorders a second array to match the first. The arrays need not be the same type.
1032676f2a66SJacob Faibussowitsch 
1033676f2a66SJacob Faibussowitsch   Not Collective
1034676f2a66SJacob Faibussowitsch 
1035676f2a66SJacob Faibussowitsch   Input Parameters:
1036676f2a66SJacob Faibussowitsch + n     - number of values
1037676f2a66SJacob Faibussowitsch . asize - size in bytes of the datatype held in arr
10386b867d5aSJose E. Roman . bsize - size in bytes of the datatype held in barr
10394d3610e3SJacob Faibussowitsch . cmp   - function pointer to comparison function
10404d3610e3SJacob Faibussowitsch - ctx   - optional context to be passed to comparison function, NULL if not needed
1041676f2a66SJacob Faibussowitsch 
10426b867d5aSJose E. Roman   Input/Output Parameters:
10436b867d5aSJose E. Roman + arr  - array to be sorted, on output it is sorted
10446b867d5aSJose E. Roman - barr - array to be reordered, on output it is reordered
10458da36234SJacob Faibussowitsch 
10468da36234SJacob Faibussowitsch   Notes:
10478da36234SJacob 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
10488da36234SJacob Faibussowitsch   overlap.
10498da36234SJacob Faibussowitsch 
10508da36234SJacob Faibussowitsch   Timsort makes the assumption that input data is already likely partially ordered, or that it contains contiguous
10518da36234SJacob Faibussowitsch   sections (termed 'runs') where the data is locally ordered (but not necessarily globally ordered). It therefore aims
10528da36234SJacob Faibussowitsch  to select slices of the array in such a way that resulting mergesorts operate on near perfectly length-balanced
10538da36234SJacob Faibussowitsch  arrays. To do so it repeatedly triggers attempts throughout to merge adjacent runs.
10548da36234SJacob Faibussowitsch 
10558da36234SJacob Faibussowitsch   Should one run continuously "win" a comparison the algorithm begins the "gallop" phase. It will aggressively search
10568da36234SJacob Faibussowitsch   the "winner" for the location of the "losers" next entry (and vice versa) to copy all preceding elements into place in
10578da36234SJacob Faibussowitsch   bulk. However if the data is truly unordered (as is the case with random data) the immense gains possible from these
10588da36234SJacob Faibussowitsch   searches are expected __not__ to repay their costs. While adjacent arrays are almost all nearly the same size, they
10598da36234SJacob Faibussowitsch   likely all contain similar data.
1060676f2a66SJacob Faibussowitsch 
1061676f2a66SJacob Faibussowitsch   Sample usage:
1062676f2a66SJacob Faibussowitsch   The comparison function must follow the qsort() comparison function paradigm, returning the sign of the difference
1063676f2a66SJacob Faibussowitsch   between its arguments. If left < right : return -1, if left == right : return 0, if left > right : return 1. The user
1064676f2a66SJacob Faibussowitsch   may also change or reverse the order of the sort by flipping the above. Note that stability of the sort is only
1065676f2a66SJacob Faibussowitsch   guaranteed if the comparison function forms a valid trigraph. For example when sorting an array of type "my_type" in
10668da36234SJacob Faibussowitsch   increasing order
1067676f2a66SJacob Faibussowitsch 
1068676f2a66SJacob Faibussowitsch .vb
10694d3610e3SJacob Faibussowitsch   int my_increasing_comparison_function(const void *left, const void *right, void *ctx) {
1070676f2a66SJacob Faibussowitsch     my_type l = *(my_type *) left, r = *(my_type *) right;
1071676f2a66SJacob Faibussowitsch     return (l < r) ? -1 : (l > r);
1072676f2a66SJacob Faibussowitsch   }
1073676f2a66SJacob Faibussowitsch .ve
10744d3610e3SJacob Faibussowitsch   Note the context is unused here but you may use it to pass and subsequently access whatever information required
10754d3610e3SJacob Faibussowitsch   inside the comparison function. The context pointer will unaltered except for any changes made inside the comparison function.
1076676f2a66SJacob Faibussowitsch   Then pass the function
1077676f2a66SJacob Faibussowitsch .vb
10784d3610e3SJacob Faibussowitsch   PetscTimSortWithArray(n, arr, sizeof(arr[0]), barr, sizeof(barr[0]), my_increasing_comparison_function, ctx)
1079676f2a66SJacob Faibussowitsch .ve
1080676f2a66SJacob Faibussowitsch 
10814d3610e3SJacob Faibussowitsch   Fortran Notes:
10824d3610e3SJacob Faibussowitsch   To use this from fortran you must write a comparison subroutine with 4 arguments which accepts left, right, ctx and
10838da36234SJacob Faibussowitsch   returns result. For example
10844d3610e3SJacob Faibussowitsch .vb
10854d3610e3SJacob Faibussowitsch  subroutine CompareIntegers(left,right,ctx,result)
10864d3610e3SJacob Faibussowitsch    implicit none
10874d3610e3SJacob Faibussowitsch 
10884d3610e3SJacob Faibussowitsch    PetscInt,intent(in) :: left, right
10894d3610e3SJacob Faibussowitsch    type(UserCtx)       :: ctx
10904d3610e3SJacob Faibussowitsch    integer,intent(out) :: result
10914d3610e3SJacob Faibussowitsch 
10924d3610e3SJacob Faibussowitsch    if (left < right) then
10934d3610e3SJacob Faibussowitsch      result = -1
10944d3610e3SJacob Faibussowitsch    else if (left == right) then
10954d3610e3SJacob Faibussowitsch      result = 0
10964d3610e3SJacob Faibussowitsch    else
10974d3610e3SJacob Faibussowitsch      result = 1
10984d3610e3SJacob Faibussowitsch    end if
10994d3610e3SJacob Faibussowitsch    return
11004d3610e3SJacob Faibussowitsch  end subroutine CompareIntegers
11014d3610e3SJacob Faibussowitsch .ve
11024d3610e3SJacob Faibussowitsch 
1103676f2a66SJacob Faibussowitsch   References:
1104676f2a66SJacob Faibussowitsch   1. - Tim Peters. https://bugs.python.org/file4451/timsort.txt
1105676f2a66SJacob Faibussowitsch 
1106676f2a66SJacob Faibussowitsch   Level: developer
1107676f2a66SJacob Faibussowitsch 
1108676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscIntSortSemiOrderedWithArray(), PetscRealSortSemiOrderedWithArrayInt(), PetscMPIIntSortSemiOrderedWithArray()
11094d3610e3SJacob Faibussowitsch @*/
11104d3610e3SJacob Faibussowitsch PetscErrorCode PetscTimSortWithArray(PetscInt n, void *arr, size_t asize, void *barr, size_t bsize, int (*cmp)(const void *, const void *, void *), void *ctx)
1111676f2a66SJacob Faibussowitsch {
1112676f2a66SJacob Faibussowitsch   PetscInt           stacksize = 0, minrun, runstart = 0, runend = 0;
1113676f2a66SJacob Faibussowitsch   PetscTimSortStack  runstack[128];
1114676f2a66SJacob Faibussowitsch   PetscTimSortBuffer abuff, bbuff;
1115676f2a66SJacob Faibussowitsch   PetscErrorCode     ierr;
1116676f2a66SJacob Faibussowitsch   /* stacksize  = log_phi(n) = log_2(n)/log_2(phi), so 128 is enough for ~5.614e26 elements.
1117676f2a66SJacob Faibussowitsch    It is so unlikely that this limit is reached that this is __never__ checked for */
1118676f2a66SJacob Faibussowitsch 
1119676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1120676f2a66SJacob Faibussowitsch   /* Compute minrun. Minrun should be (32, 65) such that N/minrun
1121676f2a66SJacob Faibussowitsch    is a power of 2 or one plus a power of 2 */
1122676f2a66SJacob Faibussowitsch   {
1123676f2a66SJacob Faibussowitsch     PetscInt t = n, r = 0;
1124676f2a66SJacob Faibussowitsch     /* r becomes 1 if the least significant bits contain at least one off bit */
1125676f2a66SJacob Faibussowitsch     while (t >= 64) {
1126676f2a66SJacob Faibussowitsch       r |= t & 1;
1127676f2a66SJacob Faibussowitsch       t >>= 1;
1128676f2a66SJacob Faibussowitsch     }
1129676f2a66SJacob Faibussowitsch     minrun = t + r;
1130676f2a66SJacob Faibussowitsch   }
11314d3610e3SJacob Faibussowitsch   if (PetscDefined(USE_DEBUG)) {
1132*3ca90d2dSJacob Faibussowitsch     ierr = PetscInfo1(NULL, "minrun = %" PetscInt_FMT "\n", minrun);CHKERRQ(ierr);
1133*3ca90d2dSJacob Faibussowitsch     if ((n >= 64) && ((minrun < 32) || (minrun > 65))) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Calculated minrun %" PetscInt_FMT " not in range (32,65)",minrun);
11344d3610e3SJacob Faibussowitsch   }
1135676f2a66SJacob Faibussowitsch   ierr = PetscMalloc1((size_t) minrun*asize, &abuff.ptr);CHKERRQ(ierr);
1136676f2a66SJacob Faibussowitsch   abuff.size = (size_t) minrun*asize;
1137676f2a66SJacob Faibussowitsch   abuff.maxsize = (size_t) n*asize;
1138676f2a66SJacob Faibussowitsch   ierr = PetscMalloc1((size_t) minrun*bsize, &bbuff.ptr);CHKERRQ(ierr);
1139676f2a66SJacob Faibussowitsch   bbuff.size = (size_t) minrun*bsize;
1140676f2a66SJacob Faibussowitsch   bbuff.maxsize = (size_t) n*bsize;
1141676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
1142676f2a66SJacob Faibussowitsch   while (runstart < n) {
1143676f2a66SJacob Faibussowitsch     /* Check if additional entries are at least partially ordered and build natural run */
11444d3610e3SJacob Faibussowitsch     ierr = PetscTimSortBuildRunWithArray_Private((char *)arr, abuff.ptr, asize, (char *)barr, bbuff.ptr, bsize, cmp, ctx, n, minrun, runstart, &runend);CHKERRQ(ierr);
1145676f2a66SJacob Faibussowitsch     runstack[stacksize].start = runstart;
1146676f2a66SJacob Faibussowitsch     runstack[stacksize].size = runend-runstart+1;
11474d3610e3SJacob Faibussowitsch     ierr = PetscTimSortMergeCollapseWithArray_Private((char *)arr, asize, (char *)barr, bsize, cmp, ctx, &abuff, &bbuff, runstack, &stacksize);CHKERRQ(ierr);
1148676f2a66SJacob Faibussowitsch     ++stacksize;
1149676f2a66SJacob Faibussowitsch     runstart = runend+1;
1150676f2a66SJacob Faibussowitsch   }
1151676f2a66SJacob Faibussowitsch   /* Have been inside while, so discard last stacksize++ */
1152676f2a66SJacob Faibussowitsch   --stacksize;
11534d3610e3SJacob Faibussowitsch   ierr = PetscTimSortForceCollapseWithArray_Private((char *)arr, asize, (char *)barr, bsize, cmp, ctx, &abuff, &bbuff, runstack, stacksize);CHKERRQ(ierr);
1154676f2a66SJacob Faibussowitsch   ierr = PetscFree(abuff.ptr);CHKERRQ(ierr);
1155676f2a66SJacob Faibussowitsch   ierr = PetscFree(bbuff.ptr);CHKERRQ(ierr);
1156676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
1157676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1158676f2a66SJacob Faibussowitsch }
1159676f2a66SJacob Faibussowitsch 
1160676f2a66SJacob Faibussowitsch /*@
1161676f2a66SJacob Faibussowitsch    PetscIntSortSemiOrdered - Sorts an array of integers in place in increasing order.
1162676f2a66SJacob Faibussowitsch 
1163676f2a66SJacob Faibussowitsch    Not Collective
1164676f2a66SJacob Faibussowitsch 
1165676f2a66SJacob Faibussowitsch    Input Parameters:
1166676f2a66SJacob Faibussowitsch +  n   - number of values
1167676f2a66SJacob Faibussowitsch -  arr - array of integers
1168676f2a66SJacob Faibussowitsch 
1169676f2a66SJacob Faibussowitsch    Output Parameters:
1170676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1171676f2a66SJacob Faibussowitsch 
1172676f2a66SJacob Faibussowitsch    Notes:
1173676f2a66SJacob Faibussowitsch    If the array is less than 64 entries long PetscSortInt() is automatically used.
1174676f2a66SJacob Faibussowitsch 
1175676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortInt(). While this function works for any array of integers it is
1176676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1177a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1178676f2a66SJacob Faibussowitsch 
1179676f2a66SJacob Faibussowitsch    Level: intermediate
1180676f2a66SJacob Faibussowitsch 
1181676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscSortInt(), PetscSortIntWithPermutation()
1182676f2a66SJacob Faibussowitsch @*/
1183676f2a66SJacob Faibussowitsch PetscErrorCode PetscIntSortSemiOrdered(PetscInt n, PetscInt arr[])
1184676f2a66SJacob Faibussowitsch {
1185676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
1186676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1187676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1188676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr,2);
1189676f2a66SJacob Faibussowitsch   if (n < 64) {
1190676f2a66SJacob Faibussowitsch     ierr = PetscSortInt(n, arr);CHKERRQ(ierr);
1191676f2a66SJacob Faibussowitsch   } else {
11924d3610e3SJacob Faibussowitsch     ierr = PetscTimSort(n, arr, sizeof(PetscInt), Compare_PetscInt_Private, NULL);CHKERRQ(ierr);
1193676f2a66SJacob Faibussowitsch   }
1194676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1195676f2a66SJacob Faibussowitsch }
1196676f2a66SJacob Faibussowitsch 
1197676f2a66SJacob Faibussowitsch /*@
1198676f2a66SJacob Faibussowitsch    PetscIntSortSemiOrderedWithArray - Sorts an array of integers in place in increasing order and reorders a second
1199676f2a66SJacob Faibussowitsch    array to match the first.
1200676f2a66SJacob Faibussowitsch 
1201676f2a66SJacob Faibussowitsch    Not Collective
1202676f2a66SJacob Faibussowitsch 
12036b867d5aSJose E. Roman    Input Parameter:
12046b867d5aSJose E. Roman .  n   - number of values
1205676f2a66SJacob Faibussowitsch 
12066b867d5aSJose E. Roman    Input/Output Parameters:
12076b867d5aSJose E. Roman +  arr1 - array of integers to be sorted, modified on output
12086b867d5aSJose E. Roman -  arr2 - array of integers to be reordered, modified on output
1209676f2a66SJacob Faibussowitsch 
1210676f2a66SJacob Faibussowitsch    Notes:
1211676f2a66SJacob Faibussowitsch    The arrays CANNOT overlap.
1212676f2a66SJacob Faibussowitsch 
1213676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortIntWithArray(). While this function works for any array of integers it is
1214676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1215a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1216676f2a66SJacob Faibussowitsch 
1217676f2a66SJacob Faibussowitsch    Level: intermediate
1218676f2a66SJacob Faibussowitsch 
1219676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscSortIntWithArray(), PetscSortIntWithPermutation()
1220676f2a66SJacob Faibussowitsch @*/
1221676f2a66SJacob Faibussowitsch PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt n, PetscInt arr1[], PetscInt arr2[])
1222676f2a66SJacob Faibussowitsch {
1223676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
1224676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
122552e73695SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1226676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr1,2);
1227676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2,3);
12286e260675SJacob Faibussowitsch   /* cannot export out to PetscIntSortWithArray here since it isn't stable */
12294d3610e3SJacob Faibussowitsch   ierr = PetscTimSortWithArray(n, arr1, sizeof(PetscInt), arr2, sizeof(PetscInt), Compare_PetscInt_Private, NULL);CHKERRQ(ierr);
1230676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1231676f2a66SJacob Faibussowitsch }
1232676f2a66SJacob Faibussowitsch 
1233676f2a66SJacob Faibussowitsch /*@
1234676f2a66SJacob Faibussowitsch    PetscMPIIntSortSemiOrdered - Sorts an array of PetscMPIInts in place in increasing order.
1235676f2a66SJacob Faibussowitsch 
1236676f2a66SJacob Faibussowitsch    Not Collective
1237676f2a66SJacob Faibussowitsch 
1238676f2a66SJacob Faibussowitsch    Input Parameters:
1239676f2a66SJacob Faibussowitsch +  n   - number of values
1240676f2a66SJacob Faibussowitsch -  arr - array of PetscMPIInts
1241676f2a66SJacob Faibussowitsch 
1242676f2a66SJacob Faibussowitsch    Output Parameters:
1243676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1244676f2a66SJacob Faibussowitsch 
1245676f2a66SJacob Faibussowitsch    Notes:
1246676f2a66SJacob Faibussowitsch    If the array is less than 64 entries long PetscSortMPIInt() is automatically used.
1247676f2a66SJacob Faibussowitsch 
1248676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortMPIInt(). While this function works for any array of PetscMPIInts it is
1249676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1250a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1251676f2a66SJacob Faibussowitsch 
1252676f2a66SJacob Faibussowitsch    Level: intermediate
1253676f2a66SJacob Faibussowitsch 
1254676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscSortMPIInt()
1255676f2a66SJacob Faibussowitsch @*/
1256676f2a66SJacob Faibussowitsch PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt n, PetscMPIInt arr[])
1257676f2a66SJacob Faibussowitsch {
1258676f2a66SJacob Faibussowitsch   PetscErrorCode  ierr;
1259676f2a66SJacob Faibussowitsch 
1260676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
126152e73695SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1262676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr,2);
1263676f2a66SJacob Faibussowitsch   if (n < 64) {
1264676f2a66SJacob Faibussowitsch     ierr = PetscSortMPIInt(n, arr);CHKERRQ(ierr);
1265676f2a66SJacob Faibussowitsch   } else {
12664d3610e3SJacob Faibussowitsch     ierr = PetscTimSort(n, arr, sizeof(PetscMPIInt), Compare_PetscMPIInt_Private, NULL);CHKERRQ(ierr);
1267676f2a66SJacob Faibussowitsch   }
1268676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1269676f2a66SJacob Faibussowitsch }
1270676f2a66SJacob Faibussowitsch 
1271676f2a66SJacob Faibussowitsch /*@
1272676f2a66SJacob Faibussowitsch    PetscMPIIntSortSemiOrderedWithArray - Sorts an array of integers in place in increasing order and reorders a second
1273676f2a66SJacob Faibussowitsch    array to match the first.
1274676f2a66SJacob Faibussowitsch 
1275676f2a66SJacob Faibussowitsch    Not Collective
1276676f2a66SJacob Faibussowitsch 
12776b867d5aSJose E. Roman    Input Parameter:
12786b867d5aSJose E. Roman .  n   - number of values
1279676f2a66SJacob Faibussowitsch 
12806b867d5aSJose E. Roman    Input/Output Parameters:
12816b867d5aSJose E. Roman .  arr1 - array of integers to be sorted, modified on output
12826b867d5aSJose E. Roman -  arr2 - array of integers to be reordered, modified on output
1283676f2a66SJacob Faibussowitsch 
1284676f2a66SJacob Faibussowitsch    Notes:
1285676f2a66SJacob Faibussowitsch    The arrays CANNOT overlap.
1286676f2a66SJacob Faibussowitsch 
1287676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortMPIIntWithArray(). While this function works for any array of integers it is
1288676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1289a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1290676f2a66SJacob Faibussowitsch 
1291676f2a66SJacob Faibussowitsch    Level: intermediate
1292676f2a66SJacob Faibussowitsch 
1293676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscSortMPIIntWithArray(), PetscSortMPIIntWithPermutation()
1294676f2a66SJacob Faibussowitsch @*/
1295676f2a66SJacob Faibussowitsch PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt n, PetscMPIInt arr1[], PetscMPIInt arr2[])
1296676f2a66SJacob Faibussowitsch {
1297676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
1298676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1299676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1300676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr1,2);
1301676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2,3);
13026e260675SJacob Faibussowitsch   /* cannot export out to PetscMPIIntSortWithArray here since it isn't stable */
13034d3610e3SJacob Faibussowitsch   ierr = PetscTimSortWithArray(n, arr1, sizeof(PetscMPIInt), arr2, sizeof(PetscMPIInt), Compare_PetscMPIInt_Private, NULL);CHKERRQ(ierr);
1304676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1305676f2a66SJacob Faibussowitsch }
1306676f2a66SJacob Faibussowitsch 
1307676f2a66SJacob Faibussowitsch /*@
1308676f2a66SJacob Faibussowitsch    PetscRealSortSemiOrdered - Sorts an array of PetscReals in place in increasing order.
1309676f2a66SJacob Faibussowitsch 
1310676f2a66SJacob Faibussowitsch    Not Collective
1311676f2a66SJacob Faibussowitsch 
1312676f2a66SJacob Faibussowitsch    Input Parameters:
1313676f2a66SJacob Faibussowitsch +  n   - number of values
1314676f2a66SJacob Faibussowitsch -  arr - array of PetscReals
1315676f2a66SJacob Faibussowitsch 
1316676f2a66SJacob Faibussowitsch    Output Parameters:
1317676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1318676f2a66SJacob Faibussowitsch 
1319676f2a66SJacob Faibussowitsch    Notes:
1320676f2a66SJacob Faibussowitsch    If the array is less than 64 entries long PetscSortReal() is automatically used.
1321676f2a66SJacob Faibussowitsch 
1322676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortReal(). While this function works for any array of PetscReals it is
1323676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1324a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1325676f2a66SJacob Faibussowitsch 
1326676f2a66SJacob Faibussowitsch    Level: intermediate
1327676f2a66SJacob Faibussowitsch 
1328676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscSortReal(), PetscSortRealWithPermutation()
1329676f2a66SJacob Faibussowitsch @*/
1330676f2a66SJacob Faibussowitsch PetscErrorCode PetscRealSortSemiOrdered(PetscInt n, PetscReal arr[])
1331676f2a66SJacob Faibussowitsch {
1332676f2a66SJacob Faibussowitsch   PetscErrorCode  ierr;
1333676f2a66SJacob Faibussowitsch 
1334676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1335676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1336676f2a66SJacob Faibussowitsch   PetscValidRealPointer(arr,2);
1337676f2a66SJacob Faibussowitsch   if (n < 64) {
1338676f2a66SJacob Faibussowitsch     ierr = PetscSortReal(n, arr);CHKERRQ(ierr);
1339676f2a66SJacob Faibussowitsch   } else {
13404d3610e3SJacob Faibussowitsch     ierr = PetscTimSort(n, arr, sizeof(PetscReal), Compare_PetscReal_Private, NULL);CHKERRQ(ierr);
1341676f2a66SJacob Faibussowitsch   }
1342676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1343676f2a66SJacob Faibussowitsch }
1344676f2a66SJacob Faibussowitsch 
1345676f2a66SJacob Faibussowitsch /*@
1346676f2a66SJacob Faibussowitsch    PetscRealSortSemiOrderedWithArrayInt - Sorts an array of PetscReals in place in increasing order and reorders a second
1347676f2a66SJacob Faibussowitsch    array of PetscInts to match the first.
1348676f2a66SJacob Faibussowitsch 
1349676f2a66SJacob Faibussowitsch    Not Collective
1350676f2a66SJacob Faibussowitsch 
13516b867d5aSJose E. Roman    Input Parameter:
13526b867d5aSJose E. Roman .  n   - number of values
1353676f2a66SJacob Faibussowitsch 
13546b867d5aSJose E. Roman    Input/Output Parameters:
13556b867d5aSJose E. Roman .  arr1 - array of PetscReals to be sorted, modified on output
13566b867d5aSJose E. Roman -  arr2 - array of PetscReals to be reordered, modified on output
1357676f2a66SJacob Faibussowitsch 
1358676f2a66SJacob Faibussowitsch    Notes:
1359676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortRealWithArray(). While this function works for any array of PetscReals it is
1360676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1361a5b23f4aSJose E. Roman    recommended that the user benchmark their code to see which routine is fastest.
1362676f2a66SJacob Faibussowitsch 
1363676f2a66SJacob Faibussowitsch    Level: intermediate
1364676f2a66SJacob Faibussowitsch 
1365676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscSortRealWithArrayInt(), PetscSortRealWithPermutation()
1366676f2a66SJacob Faibussowitsch @*/
1367676f2a66SJacob Faibussowitsch PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt n, PetscReal arr1[], PetscInt arr2[])
1368676f2a66SJacob Faibussowitsch {
1369676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
1370676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1371676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1372676f2a66SJacob Faibussowitsch   PetscValidRealPointer(arr1,2);
1373676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2,3);
13746e260675SJacob Faibussowitsch   /* cannot export out to PetscRealSortWithArrayInt here since it isn't stable */
13754d3610e3SJacob Faibussowitsch   ierr = PetscTimSortWithArray(n, arr1, sizeof(PetscReal), arr2, sizeof(PetscInt), Compare_PetscReal_Private, NULL);CHKERRQ(ierr);
1376676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1377676f2a66SJacob Faibussowitsch }
1378