xref: /petsc/src/sys/utils/sortso.c (revision 676f2a665bc48dc8c31863fa9f09556a43eff7ab)
1*676f2a66SJacob Faibussowitsch #include <petscsys.h>                /*I  "petscsys.h"  I*/
2*676f2a66SJacob Faibussowitsch #include <petsc/private/petscimpl.h>
3*676f2a66SJacob Faibussowitsch 
4*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE int Compare_PetscMPIInt_Private(const void *left, const void *right)
5*676f2a66SJacob Faibussowitsch {
6*676f2a66SJacob Faibussowitsch   PetscMPIInt l = *(PetscMPIInt *) left, r = *(PetscMPIInt *) right;
7*676f2a66SJacob Faibussowitsch   return (l < r) ? -1 : (l > r);
8*676f2a66SJacob Faibussowitsch }
9*676f2a66SJacob Faibussowitsch 
10*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE int Compare_PetscInt_Private(const void *left, const void *right)
11*676f2a66SJacob Faibussowitsch {
12*676f2a66SJacob Faibussowitsch   PetscInt l = *(PetscInt *) left, r = *(PetscInt *) right;
13*676f2a66SJacob Faibussowitsch   return (l < r) ? -1 : (l > r);
14*676f2a66SJacob Faibussowitsch }
15*676f2a66SJacob Faibussowitsch 
16*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE int Compare_PetscReal_Private(const void *left, const void *right)
17*676f2a66SJacob Faibussowitsch {
18*676f2a66SJacob Faibussowitsch   PetscReal l = *(PetscReal *) left, r = *(PetscReal *) right;
19*676f2a66SJacob Faibussowitsch   return (l < r) ? -1 : (l > r);
20*676f2a66SJacob Faibussowitsch }
21*676f2a66SJacob Faibussowitsch 
22*676f2a66SJacob Faibussowitsch #define MIN_GALLOP_CONST_GLOBAL 8
23*676f2a66SJacob Faibussowitsch static PetscInt MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
24*676f2a66SJacob Faibussowitsch typedef int (*CompFunc)(const void *, const void *);
25*676f2a66SJacob Faibussowitsch 
26*676f2a66SJacob Faibussowitsch /* Mostly to force clang uses the builtins, but can't hurt to have gcc compilers also do the same */
27*676f2a66SJacob Faibussowitsch #if !defined (PETSC_USE_DEBUG) && defined(__GNUC__)
28*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void COPYSWAPPY(char *a, char *b, char *t, size_t size)
29*676f2a66SJacob Faibussowitsch {
30*676f2a66SJacob Faibussowitsch   __builtin_memcpy(t, b, size);
31*676f2a66SJacob Faibussowitsch   __builtin_memmove(b, a, size);
32*676f2a66SJacob Faibussowitsch   __builtin_memcpy(a, t, size);
33*676f2a66SJacob Faibussowitsch   return;
34*676f2a66SJacob Faibussowitsch }
35*676f2a66SJacob Faibussowitsch 
36*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void COPYSWAPPY2(char *aa, char *ab, size_t asize, char *ba, char *bb, size_t bsize, char *t)
37*676f2a66SJacob Faibussowitsch {
38*676f2a66SJacob Faibussowitsch   __builtin_memcpy(t, ab, asize);
39*676f2a66SJacob Faibussowitsch   __builtin_memmove(ab, aa, asize);
40*676f2a66SJacob Faibussowitsch   __builtin_memcpy(aa, t, asize);
41*676f2a66SJacob Faibussowitsch   __builtin_memcpy(t, bb, bsize);
42*676f2a66SJacob Faibussowitsch   __builtin_memmove(bb, ba, bsize);
43*676f2a66SJacob Faibussowitsch   __builtin_memcpy(ba, t, bsize);
44*676f2a66SJacob Faibussowitsch   return;
45*676f2a66SJacob Faibussowitsch }
46*676f2a66SJacob Faibussowitsch 
47*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memcpy(char *dest, const char *src, size_t size)
48*676f2a66SJacob Faibussowitsch {
49*676f2a66SJacob Faibussowitsch   __builtin_memcpy(dest, src, size);
50*676f2a66SJacob Faibussowitsch   return;
51*676f2a66SJacob Faibussowitsch }
52*676f2a66SJacob Faibussowitsch 
53*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memcpy2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize)
54*676f2a66SJacob Faibussowitsch {
55*676f2a66SJacob Faibussowitsch   __builtin_memcpy(adest, asrc, asize);
56*676f2a66SJacob Faibussowitsch   __builtin_memcpy(bdest, asrc, bsize);
57*676f2a66SJacob Faibussowitsch   return;
58*676f2a66SJacob Faibussowitsch }
59*676f2a66SJacob Faibussowitsch 
60*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memmove(char *dest, const char *src, size_t size)
61*676f2a66SJacob Faibussowitsch {
62*676f2a66SJacob Faibussowitsch   __builtin_memmove(dest, src, size);
63*676f2a66SJacob Faibussowitsch   return;
64*676f2a66SJacob Faibussowitsch }
65*676f2a66SJacob Faibussowitsch 
66*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memmove2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize)
67*676f2a66SJacob Faibussowitsch {
68*676f2a66SJacob Faibussowitsch   __builtin_memmove(adest, asrc, asize);
69*676f2a66SJacob Faibussowitsch   __builtin_memmove(bdest, bsrc, bsize);
70*676f2a66SJacob Faibussowitsch   return;
71*676f2a66SJacob Faibussowitsch }
72*676f2a66SJacob Faibussowitsch # else
73*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void COPYSWAPPY(char *a, char *b, char *t, size_t size)
74*676f2a66SJacob Faibussowitsch {
75*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
76*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
77*676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(t, b, size);CHKERRABORT(PETSC_COMM_SELF,ierr);
78*676f2a66SJacob Faibussowitsch   ierr = PetscMemmove(b, a, size);CHKERRABORT(PETSC_COMM_SELF,ierr);
79*676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(a, t, size);CHKERRABORT(PETSC_COMM_SELF,ierr);
80*676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
81*676f2a66SJacob Faibussowitsch }
82*676f2a66SJacob Faibussowitsch 
83*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void COPYSWAPPY2(char *aa, char *ab, size_t asize, char *ba, char *bb, size_t bsize, char *t)
84*676f2a66SJacob Faibussowitsch {
85*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
86*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
87*676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(t, ab, asize);CHKERRABORT(PETSC_COMM_SELF,ierr);
88*676f2a66SJacob Faibussowitsch   ierr = PetscMemmove(ab, aa, asize);CHKERRABORT(PETSC_COMM_SELF,ierr);
89*676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(aa, t, asize);CHKERRABORT(PETSC_COMM_SELF,ierr);
90*676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(t, bb, bsize);CHKERRABORT(PETSC_COMM_SELF,ierr);
91*676f2a66SJacob Faibussowitsch   ierr = PetscMemmove(bb, ba, bsize);CHKERRABORT(PETSC_COMM_SELF,ierr);
92*676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(ba, t, bsize);CHKERRABORT(PETSC_COMM_SELF,ierr);
93*676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
94*676f2a66SJacob Faibussowitsch }
95*676f2a66SJacob Faibussowitsch 
96*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memcpy(char *dest, const char *src, size_t size)
97*676f2a66SJacob Faibussowitsch {
98*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
99*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
100*676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(dest, src, size);CHKERRABORT(PETSC_COMM_SELF,ierr);
101*676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
102*676f2a66SJacob Faibussowitsch }
103*676f2a66SJacob Faibussowitsch 
104*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memcpy2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize)
105*676f2a66SJacob Faibussowitsch {
106*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
107*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
108*676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(adest, asrc, asize);CHKERRABORT(PETSC_COMM_SELF,ierr);
109*676f2a66SJacob Faibussowitsch   ierr = PetscMemcpy(bdest, asrc, bsize);CHKERRABORT(PETSC_COMM_SELF,ierr);
110*676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
111*676f2a66SJacob Faibussowitsch }
112*676f2a66SJacob Faibussowitsch 
113*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memmove(char *dest, const char *src, size_t size)
114*676f2a66SJacob Faibussowitsch {
115*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
116*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
117*676f2a66SJacob Faibussowitsch   ierr = PetscMemmove(dest, src, size);CHKERRABORT(PETSC_COMM_SELF,ierr);
118*676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
119*676f2a66SJacob Faibussowitsch }
120*676f2a66SJacob Faibussowitsch 
121*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE void Petsc_memmove2(char *adest, const char *asrc, size_t asize, char *bdest, const char *bsrc, size_t bsize)
122*676f2a66SJacob Faibussowitsch {
123*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
124*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
125*676f2a66SJacob Faibussowitsch   ierr = PetscMemmove(adest, asrc, asize);CHKERRABORT(PETSC_COMM_SELF,ierr);
126*676f2a66SJacob Faibussowitsch   ierr = PetscMemmove(bdest, bsrc, bsize);CHKERRABORT(PETSC_COMM_SELF,ierr);
127*676f2a66SJacob Faibussowitsch   PetscFunctionReturnVoid();
128*676f2a66SJacob Faibussowitsch }
129*676f2a66SJacob Faibussowitsch #endif
130*676f2a66SJacob Faibussowitsch 
131*676f2a66SJacob 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] >
132*676f2a66SJacob Faibussowitsch  x. Output also inclusive.
133*676f2a66SJacob Faibussowitsch 
134*676f2a66SJacob Faibussowitsch  NOTE: Both gallopsearch functions CANNOT distinguish between inserting AFTER the entire array vs inserting at entry n!! For example for an array:
135*676f2a66SJacob Faibussowitsch 
136*676f2a66SJacob Faibussowitsch    {0,1,2,3,4,5}
137*676f2a66SJacob Faibussowitsch 
138*676f2a66SJacob Faibussowitsch    when looking to insert "5" this routine will return *m = 6, but when looking to insert "6" it will ALSO return *m = 6.
139*676f2a66SJacob Faibussowitsch   */
140*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscGallopSearchLeft_Private(const char *arr, size_t size, CompFunc cmp, PetscInt l, PetscInt r, const char *x, PetscInt *m)
141*676f2a66SJacob Faibussowitsch {
142*676f2a66SJacob Faibussowitsch   PetscInt last = l, k = 1, mid, cur = l+1;
143*676f2a66SJacob Faibussowitsch 
144*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
145*676f2a66SJacob Faibussowitsch   *m = l;
146*676f2a66SJacob Faibussowitsch   if (PetscUnlikelyDebug(r < l)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"r %D < l %D in PetscGallopSearchLeft",r,l);
147*676f2a66SJacob Faibussowitsch   if ((*cmp)(x, arr+r*size) >= 0) {*m = r; PetscFunctionReturn(0);}
148*676f2a66SJacob Faibussowitsch   if ((*cmp)(x, (arr)+l*size) < 0 || PetscUnlikely(!(r-l))) PetscFunctionReturn(0);
149*676f2a66SJacob Faibussowitsch   while (PETSC_TRUE) {
150*676f2a66SJacob Faibussowitsch     if (PetscUnlikely(cur > r)) {cur = r; break;}
151*676f2a66SJacob Faibussowitsch     if ((*cmp)(x, (arr)+cur*size) < 0) break;
152*676f2a66SJacob Faibussowitsch     last = cur;
153*676f2a66SJacob Faibussowitsch     cur += (k <<= 1) + 1;
154*676f2a66SJacob Faibussowitsch     ++k;
155*676f2a66SJacob Faibussowitsch   }
156*676f2a66SJacob Faibussowitsch   /* standard binary search but take last 0 mid 0 cur 1 into account*/
157*676f2a66SJacob Faibussowitsch   while (cur > last + 1) {
158*676f2a66SJacob Faibussowitsch     mid = last + ((cur - last) >> 1);
159*676f2a66SJacob Faibussowitsch     if ((*cmp)(x, (arr)+mid*size) < 0) {
160*676f2a66SJacob Faibussowitsch       cur = mid;
161*676f2a66SJacob Faibussowitsch     } else {
162*676f2a66SJacob Faibussowitsch       last = mid;
163*676f2a66SJacob Faibussowitsch     }
164*676f2a66SJacob Faibussowitsch   }
165*676f2a66SJacob Faibussowitsch   *m = cur;
166*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
167*676f2a66SJacob Faibussowitsch }
168*676f2a66SJacob Faibussowitsch 
169*676f2a66SJacob 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]
170*676f2a66SJacob Faibussowitsch  < x. Output also inclusive */
171*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscGallopSearchRight_Private(const char *arr, size_t size, CompFunc cmp, PetscInt l, PetscInt r, const char *x, PetscInt *m)
172*676f2a66SJacob Faibussowitsch {
173*676f2a66SJacob Faibussowitsch   PetscInt last = r, k = 1, mid, cur = r-1;
174*676f2a66SJacob Faibussowitsch 
175*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
176*676f2a66SJacob Faibussowitsch   *m = r;
177*676f2a66SJacob Faibussowitsch   if (PetscUnlikelyDebug(r < l)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"r %D < l %D in PetscGallopSearchRight",r,l);
178*676f2a66SJacob Faibussowitsch   if ((*cmp)(x, arr+l*size) <= 0) {*m = l; PetscFunctionReturn(0);}
179*676f2a66SJacob Faibussowitsch   if ((*cmp)(x, (arr)+r*size) > 0 || PetscUnlikely(!(r-l))) PetscFunctionReturn(0);
180*676f2a66SJacob Faibussowitsch   while (PETSC_TRUE) {
181*676f2a66SJacob Faibussowitsch     if (PetscUnlikely(cur < l)) {cur = l; break;}
182*676f2a66SJacob Faibussowitsch     if ((*cmp)(x, (arr)+cur*size) > 0) break;
183*676f2a66SJacob Faibussowitsch     last = cur;
184*676f2a66SJacob Faibussowitsch     cur -= (k <<= 1) + 1;
185*676f2a66SJacob Faibussowitsch     ++k;
186*676f2a66SJacob Faibussowitsch   }
187*676f2a66SJacob Faibussowitsch   /* standard binary search but take last r-1 mid r-1 cur r-2 into account*/
188*676f2a66SJacob Faibussowitsch   while (last > cur + 1) {
189*676f2a66SJacob Faibussowitsch     mid = last - ((last - cur) >> 1);
190*676f2a66SJacob Faibussowitsch     if ((*cmp)(x, (arr)+mid*size) > 0) {
191*676f2a66SJacob Faibussowitsch       cur = mid;
192*676f2a66SJacob Faibussowitsch     } else {
193*676f2a66SJacob Faibussowitsch       last = mid;
194*676f2a66SJacob Faibussowitsch     }
195*676f2a66SJacob Faibussowitsch   }
196*676f2a66SJacob Faibussowitsch   *m = cur;
197*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
198*676f2a66SJacob Faibussowitsch }
199*676f2a66SJacob Faibussowitsch 
200*676f2a66SJacob Faibussowitsch /* Mergesort where size of left half <= size of right half, so mergesort is done left to right. Arr should be pointer to
201*676f2a66SJacob Faibussowitsch  complete array, left is first index of left array, mid is first index of right array, right is last index of right
202*676f2a66SJacob Faibussowitsch  array */
203*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortMergeLo_Private(char *arr, char *tarr, size_t size, CompFunc cmp, PetscInt left, PetscInt mid, PetscInt right)
204*676f2a66SJacob Faibussowitsch {
205*676f2a66SJacob Faibussowitsch   PetscInt       i = 0, j = mid, k = left, gallopleft = 0, gallopright = 0;
206*676f2a66SJacob Faibussowitsch   const PetscInt llen = mid-left;
207*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
208*676f2a66SJacob Faibussowitsch 
209*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
210*676f2a66SJacob Faibussowitsch   Petsc_memcpy(tarr, arr+(left*size), llen*size);
211*676f2a66SJacob Faibussowitsch   while ((i < llen) && (j <= right)) {
212*676f2a66SJacob Faibussowitsch     if ((*cmp)(tarr+(i*size), arr+(j*size)) < 0) {
213*676f2a66SJacob Faibussowitsch       Petsc_memcpy(arr+(k*size), tarr+(i*size), size);
214*676f2a66SJacob Faibussowitsch       ++k;
215*676f2a66SJacob Faibussowitsch       gallopright = 0;
216*676f2a66SJacob Faibussowitsch       if (++i < llen && ++gallopleft >= MIN_GALLOP_GLOBAL) {
217*676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
218*676f2a66SJacob Faibussowitsch         do {
219*676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
220*676f2a66SJacob Faibussowitsch           /* search temp for right[j], can move up to that of temp into arr immediately */
221*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(tarr, size, cmp, i, llen-1, arr+(j*size), &l1);CHKERRQ(ierr);
222*676f2a66SJacob Faibussowitsch           diff1 = l1-i;
223*676f2a66SJacob Faibussowitsch           Petsc_memcpy(arr+(k*size), tarr+(i*size), diff1*size);
224*676f2a66SJacob Faibussowitsch           k += diff1;
225*676f2a66SJacob Faibussowitsch           i = l1;
226*676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
227*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, size, cmp, j, right, tarr+(i*size), &l2);CHKERRQ(ierr);
228*676f2a66SJacob Faibussowitsch           diff2 = l2-j;
229*676f2a66SJacob Faibussowitsch           Petsc_memmove((arr)+k*size, (arr)+j*size, diff2*size);
230*676f2a66SJacob Faibussowitsch           k += diff2;
231*676f2a66SJacob Faibussowitsch           j = l2;
232*676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
233*676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
234*676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
235*676f2a66SJacob Faibussowitsch       }
236*676f2a66SJacob Faibussowitsch     } else {
237*676f2a66SJacob Faibussowitsch       Petsc_memmove(arr+(k*size), arr+(j*size), size);
238*676f2a66SJacob Faibussowitsch       ++k;
239*676f2a66SJacob Faibussowitsch       gallopleft = 0;
240*676f2a66SJacob Faibussowitsch       if (++j <= right && ++gallopright >= MIN_GALLOP_GLOBAL) {
241*676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
242*676f2a66SJacob Faibussowitsch         do {
243*676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
244*676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
245*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, size, cmp, j, right, tarr+(i*size), &l2);CHKERRQ(ierr);
246*676f2a66SJacob Faibussowitsch           diff2 = l2-j;
247*676f2a66SJacob Faibussowitsch           Petsc_memmove(arr+(k*size), arr+(j*size), diff2*size);
248*676f2a66SJacob Faibussowitsch           k += diff2;
249*676f2a66SJacob Faibussowitsch           j = l2;
250*676f2a66SJacob Faibussowitsch           /* search temp for right[j], can copy up to that of temp into arr immediately */
251*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(tarr, size, cmp, i, llen-1, arr+(j*size), &l1);CHKERRQ(ierr);
252*676f2a66SJacob Faibussowitsch           diff1 = l1-i;
253*676f2a66SJacob Faibussowitsch           Petsc_memcpy(arr+(k*size), tarr+(i*size), diff1*size);
254*676f2a66SJacob Faibussowitsch           k += diff1;
255*676f2a66SJacob Faibussowitsch           i = l1;
256*676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
257*676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
258*676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
259*676f2a66SJacob Faibussowitsch       }
260*676f2a66SJacob Faibussowitsch     }
261*676f2a66SJacob Faibussowitsch   }
262*676f2a66SJacob Faibussowitsch   if (i<llen) {Petsc_memcpy(arr+(k*size), tarr+(i*size), (llen-i)*size);}
263*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
264*676f2a66SJacob Faibussowitsch }
265*676f2a66SJacob Faibussowitsch 
266*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortMergeLoWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, PetscInt left, PetscInt mid, PetscInt right)
267*676f2a66SJacob Faibussowitsch {
268*676f2a66SJacob Faibussowitsch   PetscInt       i = 0, j = mid, k = left, gallopleft = 0, gallopright = 0;
269*676f2a66SJacob Faibussowitsch   const PetscInt llen = mid-left;
270*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
271*676f2a66SJacob Faibussowitsch 
272*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
273*676f2a66SJacob Faibussowitsch   Petsc_memcpy2(atarr, arr+(left*asize), llen*asize, btarr, barr+(left*bsize), llen*bsize);
274*676f2a66SJacob Faibussowitsch   while ((i < llen) && (j <= right)) {
275*676f2a66SJacob Faibussowitsch     if ((*cmp)(atarr+(i*asize), arr+(j*asize)) < 0) {
276*676f2a66SJacob Faibussowitsch       Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), asize, barr+(k*bsize), btarr+(i*bsize), bsize);
277*676f2a66SJacob Faibussowitsch       ++k;
278*676f2a66SJacob Faibussowitsch       gallopright = 0;
279*676f2a66SJacob Faibussowitsch       if (++i < llen && ++gallopleft >= MIN_GALLOP_GLOBAL) {
280*676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
281*676f2a66SJacob Faibussowitsch         do {
282*676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
283*676f2a66SJacob Faibussowitsch           /* search temp for right[j], can move up to that of temp into arr immediately */
284*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(atarr, asize, cmp, i, llen-1, arr+(j*asize), &l1);CHKERRQ(ierr);
285*676f2a66SJacob Faibussowitsch           diff1 = l1-i;
286*676f2a66SJacob Faibussowitsch           Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), diff1*asize, barr+(k*bsize), btarr+(i*bsize), diff1*bsize);
287*676f2a66SJacob Faibussowitsch           k += diff1;
288*676f2a66SJacob Faibussowitsch           i = l1;
289*676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
290*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, j, right, atarr+(i*asize), &l2);CHKERRQ(ierr);
291*676f2a66SJacob Faibussowitsch           diff2 = l2-j;
292*676f2a66SJacob Faibussowitsch           Petsc_memmove2(arr+(k*asize), arr+(j*asize), diff2*asize, barr+(k*bsize), barr+(j*bsize), diff2*bsize);
293*676f2a66SJacob Faibussowitsch           k += diff2;
294*676f2a66SJacob Faibussowitsch           j = l2;
295*676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
296*676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
297*676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
298*676f2a66SJacob Faibussowitsch       }
299*676f2a66SJacob Faibussowitsch     } else {
300*676f2a66SJacob Faibussowitsch       Petsc_memmove2(arr+(k*asize), arr+(j*asize), asize, barr+(k*bsize), barr+(j*bsize), bsize);
301*676f2a66SJacob Faibussowitsch       ++k;
302*676f2a66SJacob Faibussowitsch       gallopleft = 0;
303*676f2a66SJacob Faibussowitsch       if (++j <= right && ++gallopright >= MIN_GALLOP_GLOBAL) {
304*676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
305*676f2a66SJacob Faibussowitsch         do {
306*676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
307*676f2a66SJacob Faibussowitsch           /* search right for temp[i], can move up to that many of right into arr */
308*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, j, right, atarr+(i*asize), &l2);CHKERRQ(ierr);
309*676f2a66SJacob Faibussowitsch           diff2 = l2-j;
310*676f2a66SJacob Faibussowitsch           Petsc_memmove2(arr+(k*asize), arr+(j*asize), diff2*asize, barr+(k*bsize), barr+(j*bsize), diff2*bsize);
311*676f2a66SJacob Faibussowitsch           k += diff2;
312*676f2a66SJacob Faibussowitsch           j = l2;
313*676f2a66SJacob Faibussowitsch           /* search temp for right[j], can copy up to that of temp into arr immediately */
314*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(atarr, asize, cmp, i, llen-1, arr+(j*asize), &l1);CHKERRQ(ierr);
315*676f2a66SJacob Faibussowitsch           diff1 = l1-i;
316*676f2a66SJacob Faibussowitsch           Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), diff1*asize, barr+(k*bsize), btarr+(i*bsize), diff1*bsize);
317*676f2a66SJacob Faibussowitsch           k += diff1;
318*676f2a66SJacob Faibussowitsch           i = l1;
319*676f2a66SJacob Faibussowitsch           if (i >= llen || j > right) break;
320*676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
321*676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
322*676f2a66SJacob Faibussowitsch       }
323*676f2a66SJacob Faibussowitsch     }
324*676f2a66SJacob Faibussowitsch   }
325*676f2a66SJacob Faibussowitsch   if (i<llen) {Petsc_memcpy2(arr+(k*asize), atarr+(i*asize), (llen-i)*asize, barr+(k*bsize), btarr+(i*bsize), (llen-i)*bsize);}
326*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
327*676f2a66SJacob Faibussowitsch }
328*676f2a66SJacob Faibussowitsch 
329*676f2a66SJacob Faibussowitsch /* Mergesort where size of right half < size of left half, so mergesort is done right to left. Arr should be pointer to
330*676f2a66SJacob Faibussowitsch  complete array, left is first index of left array, mid is first index of right array, right is last index of right
331*676f2a66SJacob Faibussowitsch  array */
332*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortMergeHi_Private(char *arr, char *tarr, size_t size, CompFunc cmp, PetscInt left, PetscInt mid, PetscInt right)
333*676f2a66SJacob Faibussowitsch {
334*676f2a66SJacob Faibussowitsch   PetscInt       i = right-mid, j = mid-1, k = right, gallopleft = 0, gallopright = 0;
335*676f2a66SJacob Faibussowitsch   const PetscInt rlen = right-mid+1;
336*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
337*676f2a66SJacob Faibussowitsch 
338*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
339*676f2a66SJacob Faibussowitsch   Petsc_memcpy(tarr, (arr)+mid*size, rlen*size);
340*676f2a66SJacob Faibussowitsch   while ((i >= 0) && (j >= left)) {
341*676f2a66SJacob Faibussowitsch     if ((*cmp)((tarr)+i*size, (arr)+j*size) > 0) {
342*676f2a66SJacob Faibussowitsch       Petsc_memcpy((arr)+k*size, (tarr)+i*size, size);
343*676f2a66SJacob Faibussowitsch       --k;
344*676f2a66SJacob Faibussowitsch       gallopleft = 0;
345*676f2a66SJacob Faibussowitsch       if (--i >= 0 && ++gallopright >= MIN_GALLOP_GLOBAL) {
346*676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
347*676f2a66SJacob Faibussowitsch         do {
348*676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
349*676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
350*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(tarr, size, cmp, 0, i, (arr)+j*size, &l1);CHKERRQ(ierr);
351*676f2a66SJacob Faibussowitsch           diff1 = i-l1;
352*676f2a66SJacob Faibussowitsch           Petsc_memcpy((arr)+(k-diff1+1)*size, (tarr)+(l1+1)*size, diff1*size);
353*676f2a66SJacob Faibussowitsch           k -= diff1;
354*676f2a66SJacob Faibussowitsch           i = l1;
355*676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
356*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, size, cmp, left, j, (tarr)+i*size, &l2);CHKERRQ(ierr);
357*676f2a66SJacob Faibussowitsch           diff2 = j-l2;
358*676f2a66SJacob Faibussowitsch           Petsc_memmove((arr)+(k-diff2+1)*size, (arr)+(l2+1)*size, diff2*size);
359*676f2a66SJacob Faibussowitsch           k -= diff2;
360*676f2a66SJacob Faibussowitsch           j = l2;
361*676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
362*676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
363*676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
364*676f2a66SJacob Faibussowitsch       }
365*676f2a66SJacob Faibussowitsch     } else {
366*676f2a66SJacob Faibussowitsch       Petsc_memmove((arr)+k*size, (arr)+j*size, size);
367*676f2a66SJacob Faibussowitsch       --k;
368*676f2a66SJacob Faibussowitsch       gallopright = 0;
369*676f2a66SJacob Faibussowitsch       if (--j >= left && ++gallopleft >= MIN_GALLOP_GLOBAL) {
370*676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
371*676f2a66SJacob Faibussowitsch         do {
372*676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
373*676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
374*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, size, cmp, left, j, (tarr)+i*size, &l2);CHKERRQ(ierr);
375*676f2a66SJacob Faibussowitsch           diff2 = j-l2;
376*676f2a66SJacob Faibussowitsch           Petsc_memmove((arr)+(k-diff2+1)*size, (arr)+(l2+1)*size, diff2*size);
377*676f2a66SJacob Faibussowitsch           k -= diff2;
378*676f2a66SJacob Faibussowitsch           j = l2;
379*676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
380*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(tarr, size, cmp, 0, i, (arr)+j*size, &l1);CHKERRQ(ierr);
381*676f2a66SJacob Faibussowitsch           diff1 = i-l1;
382*676f2a66SJacob Faibussowitsch           Petsc_memcpy((arr)+(k-diff1+1)*size, (tarr)+(l1+1)*size, diff1*size);
383*676f2a66SJacob Faibussowitsch           k -= diff1;
384*676f2a66SJacob Faibussowitsch           i = l1;
385*676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
386*676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
387*676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
388*676f2a66SJacob Faibussowitsch       }
389*676f2a66SJacob Faibussowitsch     }
390*676f2a66SJacob Faibussowitsch   }
391*676f2a66SJacob Faibussowitsch   if (i >= 0) {Petsc_memcpy((arr)+left*size, tarr, (i+1)*size);}
392*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
393*676f2a66SJacob Faibussowitsch }
394*676f2a66SJacob Faibussowitsch 
395*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortMergeHiWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, PetscInt left, PetscInt mid, PetscInt right)
396*676f2a66SJacob Faibussowitsch {
397*676f2a66SJacob Faibussowitsch   PetscInt       i = right-mid, j = mid-1, k = right, gallopleft = 0, gallopright = 0;
398*676f2a66SJacob Faibussowitsch   const PetscInt rlen = right-mid+1;
399*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
400*676f2a66SJacob Faibussowitsch 
401*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
402*676f2a66SJacob Faibussowitsch   Petsc_memcpy2(atarr, (arr)+mid*asize, rlen*asize, btarr, (barr)+mid*bsize, rlen*bsize);
403*676f2a66SJacob Faibussowitsch   while ((i >= 0) && (j >= left)) {
404*676f2a66SJacob Faibussowitsch     if ((*cmp)((atarr)+i*asize, (arr)+j*asize) > 0) {
405*676f2a66SJacob Faibussowitsch       Petsc_memcpy2((arr)+k*asize, (atarr)+i*asize, asize, (barr)+k*bsize, (btarr)+i*bsize, bsize);
406*676f2a66SJacob Faibussowitsch       --k;
407*676f2a66SJacob Faibussowitsch       gallopleft = 0;
408*676f2a66SJacob Faibussowitsch       if (--i >= 0 && ++gallopright >= MIN_GALLOP_GLOBAL) {
409*676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
410*676f2a66SJacob Faibussowitsch         do {
411*676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
412*676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
413*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(atarr, asize, cmp, 0, i, (arr)+j*asize, &l1);CHKERRQ(ierr);
414*676f2a66SJacob Faibussowitsch           diff1 = i-l1;
415*676f2a66SJacob 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);
416*676f2a66SJacob Faibussowitsch           k -= diff1;
417*676f2a66SJacob Faibussowitsch           i = l1;
418*676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
419*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, asize, cmp, left, j, (atarr)+i*asize, &l2);CHKERRQ(ierr);
420*676f2a66SJacob Faibussowitsch           diff2 = j-l2;
421*676f2a66SJacob 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);
422*676f2a66SJacob Faibussowitsch           k -= diff2;
423*676f2a66SJacob Faibussowitsch           j = l2;
424*676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
425*676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
426*676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
427*676f2a66SJacob Faibussowitsch       }
428*676f2a66SJacob Faibussowitsch     } else {
429*676f2a66SJacob Faibussowitsch       Petsc_memmove2((arr)+k*asize, (arr)+j*asize, asize, (barr)+k*bsize, (barr)+j*bsize, bsize);
430*676f2a66SJacob Faibussowitsch       --k;
431*676f2a66SJacob Faibussowitsch       gallopright = 0;
432*676f2a66SJacob Faibussowitsch       if (--j >= left && ++gallopleft >= MIN_GALLOP_GLOBAL) {
433*676f2a66SJacob Faibussowitsch         PetscInt l1, l2, diff1, diff2;
434*676f2a66SJacob Faibussowitsch         do {
435*676f2a66SJacob Faibussowitsch           if (MIN_GALLOP_GLOBAL > 1) --MIN_GALLOP_GLOBAL;
436*676f2a66SJacob Faibussowitsch           /* search left for temp[i], can move up to that many of left up arr */
437*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, asize, cmp, left, j, (atarr)+i*asize, &l2);CHKERRQ(ierr);
438*676f2a66SJacob Faibussowitsch           diff2 = j-l2;
439*676f2a66SJacob 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);
440*676f2a66SJacob Faibussowitsch           k -= diff2;
441*676f2a66SJacob Faibussowitsch           j = l2;
442*676f2a66SJacob Faibussowitsch           /* search temp for left[j], can copy up to that many of temp into arr */
443*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(atarr, asize, cmp, 0, i, (arr)+j*asize, &l1);CHKERRQ(ierr);
444*676f2a66SJacob Faibussowitsch           diff1 = i-l1;
445*676f2a66SJacob 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);
446*676f2a66SJacob Faibussowitsch           k -= diff1;
447*676f2a66SJacob Faibussowitsch           i = l1;
448*676f2a66SJacob Faibussowitsch           if (i < 0 || j < left) break;
449*676f2a66SJacob Faibussowitsch         } while (diff1 > MIN_GALLOP_GLOBAL || diff2 > MIN_GALLOP_GLOBAL);
450*676f2a66SJacob Faibussowitsch         ++MIN_GALLOP_GLOBAL;
451*676f2a66SJacob Faibussowitsch       }
452*676f2a66SJacob Faibussowitsch     }
453*676f2a66SJacob Faibussowitsch   }
454*676f2a66SJacob Faibussowitsch   if (i >= 0) {Petsc_memcpy2((arr)+left*asize, atarr, (i+1)*asize, (barr)+left*bsize, btarr, (i+1)*bsize);}
455*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
456*676f2a66SJacob Faibussowitsch }
457*676f2a66SJacob Faibussowitsch 
458*676f2a66SJacob Faibussowitsch /* Left is inclusive lower bound of array slice, start is start location of unsorted section, right is inclusive upper
459*676f2a66SJacob Faibussowitsch  bound of array slice. If unsure of where unsorted section starts or if entire length is unsorted pass start = left */
460*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscInsertionSort_Private(char *arr, char *tarr, size_t size, CompFunc cmp, PetscInt left, PetscInt start, PetscInt right)
461*676f2a66SJacob Faibussowitsch {
462*676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
463*676f2a66SJacob Faibussowitsch 
464*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
465*676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
466*676f2a66SJacob Faibussowitsch     PetscInt j = i-1;
467*676f2a66SJacob Faibussowitsch     Petsc_memcpy(tarr, arr+(i*size), size);
468*676f2a66SJacob Faibussowitsch     while ((j >= left) && ((*cmp)(tarr, (arr)+j*size) < 0)) {
469*676f2a66SJacob Faibussowitsch       COPYSWAPPY(arr+(j+1)*size, arr+j*size, tarr+size, size);
470*676f2a66SJacob Faibussowitsch       --j;
471*676f2a66SJacob Faibussowitsch     }
472*676f2a66SJacob Faibussowitsch     Petsc_memcpy((arr)+(j+1)*size, tarr, size);
473*676f2a66SJacob Faibussowitsch   }
474*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
475*676f2a66SJacob Faibussowitsch }
476*676f2a66SJacob Faibussowitsch 
477*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscInsertionSortWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, PetscInt left, PetscInt start, PetscInt right)
478*676f2a66SJacob Faibussowitsch {
479*676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
480*676f2a66SJacob Faibussowitsch 
481*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
482*676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
483*676f2a66SJacob Faibussowitsch     PetscInt j = i-1;
484*676f2a66SJacob Faibussowitsch     Petsc_memcpy2(atarr, arr+(i*asize), asize, btarr, barr+(i*bsize), bsize);
485*676f2a66SJacob Faibussowitsch     while ((j >= left) && ((*cmp)(atarr, arr+(j*asize)) < 0)) {
486*676f2a66SJacob Faibussowitsch       COPYSWAPPY2(arr+(j+1)*asize, arr+j*asize, asize, barr+(j+1)*bsize, barr+j*bsize, bsize, atarr+asize);
487*676f2a66SJacob Faibussowitsch       --j;
488*676f2a66SJacob Faibussowitsch     }
489*676f2a66SJacob Faibussowitsch     Petsc_memcpy2(arr+(j+1)*asize, atarr, asize, barr+(j+1)*bsize, btarr, bsize);
490*676f2a66SJacob Faibussowitsch   }
491*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
492*676f2a66SJacob Faibussowitsch }
493*676f2a66SJacob Faibussowitsch 
494*676f2a66SJacob Faibussowitsch /* See PetscInsertionSort_Private */
495*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscBinaryInsertionSort_Private(char *arr, char *tarr, size_t size, CompFunc cmp, PetscInt left, PetscInt start, PetscInt right)
496*676f2a66SJacob Faibussowitsch {
497*676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
498*676f2a66SJacob Faibussowitsch 
499*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
500*676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
501*676f2a66SJacob Faibussowitsch     PetscInt l = left, r = i;
502*676f2a66SJacob Faibussowitsch     Petsc_memcpy(tarr, arr+(i*size), size);
503*676f2a66SJacob Faibussowitsch     do {
504*676f2a66SJacob Faibussowitsch       const PetscInt m = l + ((r - l) >> 1);
505*676f2a66SJacob Faibussowitsch       if ((*cmp)(tarr, arr+(m*size)) < 0) {
506*676f2a66SJacob Faibussowitsch         r = m;
507*676f2a66SJacob Faibussowitsch       } else {
508*676f2a66SJacob Faibussowitsch         l = m + 1;
509*676f2a66SJacob Faibussowitsch       }
510*676f2a66SJacob Faibussowitsch     } while (l < r);
511*676f2a66SJacob Faibussowitsch     Petsc_memmove(arr+((l+1)*size), arr+(l*size), (i-l)*size);
512*676f2a66SJacob Faibussowitsch     Petsc_memcpy(arr+(l*size), tarr, size);
513*676f2a66SJacob Faibussowitsch   }
514*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
515*676f2a66SJacob Faibussowitsch }
516*676f2a66SJacob Faibussowitsch 
517*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscBinaryInsertionSortWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, PetscInt left, PetscInt start, PetscInt right)
518*676f2a66SJacob Faibussowitsch {
519*676f2a66SJacob Faibussowitsch   PetscInt i = start == left ? start+1 : start;
520*676f2a66SJacob Faibussowitsch 
521*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
522*676f2a66SJacob Faibussowitsch   for (; i <= right; ++i) {
523*676f2a66SJacob Faibussowitsch     PetscInt l = left, r = i;
524*676f2a66SJacob Faibussowitsch     Petsc_memcpy2(atarr, arr+(i*asize), asize, btarr, barr+(i*bsize), bsize);
525*676f2a66SJacob Faibussowitsch     do {
526*676f2a66SJacob Faibussowitsch       const PetscInt m = l + ((r - l) >> 1);
527*676f2a66SJacob Faibussowitsch       if ((*cmp)(atarr, arr+(m*asize)) < 0) {
528*676f2a66SJacob Faibussowitsch         r = m;
529*676f2a66SJacob Faibussowitsch       } else {
530*676f2a66SJacob Faibussowitsch         l = m + 1;
531*676f2a66SJacob Faibussowitsch       }
532*676f2a66SJacob Faibussowitsch     } while (l < r);
533*676f2a66SJacob Faibussowitsch     Petsc_memmove2(arr+((l+1)*asize), arr+(l*asize), (i-l)*asize, barr+((l+1)*bsize), barr+(l*bsize), (i-l)*bsize);
534*676f2a66SJacob Faibussowitsch     Petsc_memcpy2(arr+(l*asize), atarr, asize, barr+(l*bsize), btarr, bsize);
535*676f2a66SJacob Faibussowitsch   }
536*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
537*676f2a66SJacob Faibussowitsch }
538*676f2a66SJacob Faibussowitsch 
539*676f2a66SJacob Faibussowitsch typedef struct {
540*676f2a66SJacob Faibussowitsch   PetscInt size;
541*676f2a66SJacob Faibussowitsch   PetscInt start;
542*676f2a66SJacob Faibussowitsch } PetscTimSortStack PETSC_ATTRIBUTEALIGNED(2*sizeof(PetscInt));
543*676f2a66SJacob Faibussowitsch 
544*676f2a66SJacob Faibussowitsch typedef struct {
545*676f2a66SJacob Faibussowitsch   char   *ptr PETSC_ATTRIBUTEALIGNED(PETSC_MEMALIGN);
546*676f2a66SJacob Faibussowitsch   size_t size;
547*676f2a66SJacob Faibussowitsch   size_t maxsize;
548*676f2a66SJacob Faibussowitsch } PetscTimSortBuffer;
549*676f2a66SJacob Faibussowitsch 
550*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortResizeBuffer_Private(PetscTimSortBuffer *buff, size_t newSize)
551*676f2a66SJacob Faibussowitsch {
552*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
553*676f2a66SJacob Faibussowitsch   if (PetscLikely(newSize <= buff->size)) PetscFunctionReturn(0);
554*676f2a66SJacob Faibussowitsch   {
555*676f2a66SJacob Faibussowitsch     /* Can't be larger than n, there is merit to simply allocating buff to n to begin with */
556*676f2a66SJacob Faibussowitsch     PetscErrorCode ierr;
557*676f2a66SJacob Faibussowitsch     size_t         newMax = PetscMin(newSize*newSize, buff->maxsize);
558*676f2a66SJacob Faibussowitsch     ierr = PetscFree(buff->ptr);CHKERRQ(ierr);
559*676f2a66SJacob Faibussowitsch     ierr = PetscMalloc1(newMax, &buff->ptr);CHKERRQ(ierr);
560*676f2a66SJacob Faibussowitsch     buff->size = newMax;
561*676f2a66SJacob Faibussowitsch   }
562*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
563*676f2a66SJacob Faibussowitsch }
564*676f2a66SJacob Faibussowitsch 
565*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortForceCollapse_Private(char *arr, size_t size, CompFunc cmp, PetscTimSortBuffer *buff, PetscTimSortStack *stack, PetscInt stacksize)
566*676f2a66SJacob Faibussowitsch {
567*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
568*676f2a66SJacob Faibussowitsch   for (;stacksize; --stacksize) {
569*676f2a66SJacob Faibussowitsch     /* A = stack[i-1], B = stack[i], if A[-1] <= B[0] means sorted */
570*676f2a66SJacob Faibussowitsch     if ((*cmp)(arr+(stack[stacksize].start-1)*size, arr+(stack[stacksize].start)*size) > 0) {
571*676f2a66SJacob Faibussowitsch       PetscInt       l, m = stack[stacksize].start, r;
572*676f2a66SJacob Faibussowitsch       PetscErrorCode ierr;
573*676f2a66SJacob Faibussowitsch       /* Search A for B[0] insertion */
574*676f2a66SJacob Faibussowitsch       ierr = PetscGallopSearchLeft_Private(arr, size, cmp, stack[stacksize-1].start, stack[stacksize].start-1, (arr)+(stack[stacksize].start)*size, &l);CHKERRQ(ierr);
575*676f2a66SJacob Faibussowitsch       /* Search B for A[-1] insertion */
576*676f2a66SJacob Faibussowitsch       ierr = PetscGallopSearchRight_Private(arr, size, cmp, stack[stacksize].start, stack[stacksize].start+stack[stacksize].size-1, (arr)+(stack[stacksize].start-1)*size, &r);CHKERRQ(ierr);
577*676f2a66SJacob Faibussowitsch       if (m-l <= r-m) {
578*676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size);CHKERRQ(ierr);
579*676f2a66SJacob Faibussowitsch         ierr = PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, l, m, r);CHKERRQ(ierr);
580*676f2a66SJacob Faibussowitsch       } else {
581*676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size);CHKERRQ(ierr);
582*676f2a66SJacob Faibussowitsch         ierr = PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, l, m, r);CHKERRQ(ierr);
583*676f2a66SJacob Faibussowitsch       }
584*676f2a66SJacob Faibussowitsch     }
585*676f2a66SJacob Faibussowitsch     /* Update A with merge */
586*676f2a66SJacob Faibussowitsch     stack[stacksize-1].size += stack[stacksize].size;
587*676f2a66SJacob Faibussowitsch   }
588*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
589*676f2a66SJacob Faibussowitsch }
590*676f2a66SJacob Faibussowitsch 
591*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortForceCollapseWithArray_Private(char *arr, size_t asize, char *barr, size_t bsize, CompFunc cmp, PetscTimSortBuffer *abuff, PetscTimSortBuffer *bbuff, PetscTimSortStack *stack, PetscInt stacksize)
592*676f2a66SJacob Faibussowitsch {
593*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
594*676f2a66SJacob Faibussowitsch   for (;stacksize; --stacksize) {
595*676f2a66SJacob Faibussowitsch     /* A = stack[i-1], B = stack[i], if A[-1] <= B[0] means sorted */
596*676f2a66SJacob Faibussowitsch     if ((*cmp)(arr+(stack[stacksize].start-1)*asize, arr+(stack[stacksize].start)*asize) > 0) {
597*676f2a66SJacob Faibussowitsch       PetscInt       l, m = stack[stacksize].start, r;
598*676f2a66SJacob Faibussowitsch       PetscErrorCode ierr;
599*676f2a66SJacob Faibussowitsch       /* Search A for B[0] insertion */
600*676f2a66SJacob Faibussowitsch       ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, stack[stacksize-1].start, stack[stacksize].start-1, (arr)+(stack[stacksize].start)*asize, &l);CHKERRQ(ierr);
601*676f2a66SJacob Faibussowitsch       /* Search B for A[-1] insertion */
602*676f2a66SJacob Faibussowitsch       ierr = PetscGallopSearchRight_Private(arr, asize, cmp, stack[stacksize].start, stack[stacksize].start+stack[stacksize].size-1, (arr)+(stack[stacksize].start-1)*asize, &r);CHKERRQ(ierr);
603*676f2a66SJacob Faibussowitsch       if (m-l <= r-m) {
604*676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize);CHKERRQ(ierr);
605*676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize);CHKERRQ(ierr);
606*676f2a66SJacob Faibussowitsch         ierr = PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, l, m, r);CHKERRQ(ierr);
607*676f2a66SJacob Faibussowitsch       } else {
608*676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize);CHKERRQ(ierr);
609*676f2a66SJacob Faibussowitsch         ierr = PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize);CHKERRQ(ierr);
610*676f2a66SJacob Faibussowitsch         ierr = PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, l, m, r);CHKERRQ(ierr);
611*676f2a66SJacob Faibussowitsch       }
612*676f2a66SJacob Faibussowitsch     }
613*676f2a66SJacob Faibussowitsch     /* Update A with merge */
614*676f2a66SJacob Faibussowitsch     stack[stacksize-1].size += stack[stacksize].size;
615*676f2a66SJacob Faibussowitsch   }
616*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
617*676f2a66SJacob Faibussowitsch }
618*676f2a66SJacob Faibussowitsch 
619*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortMergeCollapse_Private(char *arr, size_t size, CompFunc cmp, PetscTimSortBuffer *buff, PetscTimSortStack *stack, PetscInt *stacksize)
620*676f2a66SJacob Faibussowitsch {
621*676f2a66SJacob Faibussowitsch   PetscInt       i = *stacksize;
622*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
623*676f2a66SJacob Faibussowitsch 
624*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
625*676f2a66SJacob Faibussowitsch   while (i) {
626*676f2a66SJacob Faibussowitsch     PetscInt l, m, r, itemp = i;
627*676f2a66SJacob Faibussowitsch 
628*676f2a66SJacob Faibussowitsch     if (i == 1) {
629*676f2a66SJacob Faibussowitsch       /* A = stack[i-1], B = stack[i] */
630*676f2a66SJacob Faibussowitsch       if (stack[i-1].size < stack[i].size) {
631*676f2a66SJacob Faibussowitsch         /* if A[-1] <= B[0] then sorted */
632*676f2a66SJacob Faibussowitsch         if ((*cmp)(arr+(stack[i].start-1)*size, arr+(stack[i].start)*size) > 0) {
633*676f2a66SJacob Faibussowitsch           m = stack[i].start;
634*676f2a66SJacob Faibussowitsch           /* Search A for B[0] insertion */
635*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, size, cmp, stack[i-1].start, stack[i].start-1, (arr)+stack[i].start*size, &l);CHKERRQ(ierr);
636*676f2a66SJacob Faibussowitsch           /* Search B for A[-1] insertion */
637*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, size, cmp, stack[i].start, stack[i].start+stack[i].size-1, arr+(stack[i].start-1)*size, &r);CHKERRQ(ierr);
638*676f2a66SJacob Faibussowitsch           if (m-l <= r-m) {
639*676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size);CHKERRQ(ierr);
640*676f2a66SJacob Faibussowitsch             ierr = PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, l, m, r);CHKERRQ(ierr);
641*676f2a66SJacob Faibussowitsch           } else {
642*676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size);CHKERRQ(ierr);
643*676f2a66SJacob Faibussowitsch             ierr = PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, l, m, r);CHKERRQ(ierr);
644*676f2a66SJacob Faibussowitsch           }
645*676f2a66SJacob Faibussowitsch         }
646*676f2a66SJacob Faibussowitsch         /* Update A with merge */
647*676f2a66SJacob Faibussowitsch         stack[i-1].size += stack[i].size;
648*676f2a66SJacob Faibussowitsch         --i;
649*676f2a66SJacob Faibussowitsch       }
650*676f2a66SJacob Faibussowitsch     } else {
651*676f2a66SJacob Faibussowitsch       /* i > 2, i.e. C exists
652*676f2a66SJacob Faibussowitsch        A = stack[i-2], B = stack[i-1], C = stack[i]; */
653*676f2a66SJacob Faibussowitsch       if (stack[i-2].size <= stack[i-1].size+stack[i].size) {
654*676f2a66SJacob Faibussowitsch         if (stack[i-2].size < stack[i].size) {
655*676f2a66SJacob Faibussowitsch           /* merge B into A, but if A[-1] <= B[0] then already sorted */
656*676f2a66SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i-1].start-1)*size, arr+(stack[i-1].start)*size) > 0) {
657*676f2a66SJacob Faibussowitsch             m = stack[i-1].start;
658*676f2a66SJacob Faibussowitsch             /* Search A for B[0] insertion */
659*676f2a66SJacob Faibussowitsch             ierr = PetscGallopSearchLeft_Private(arr, size, cmp, stack[i-2].start, stack[i-1].start-1, (arr)+(stack[i-1].start)*size, &l);CHKERRQ(ierr);
660*676f2a66SJacob Faibussowitsch             /* Search B for A[-1] insertion */
661*676f2a66SJacob Faibussowitsch             ierr = PetscGallopSearchRight_Private(arr, size, cmp, stack[i-1].start, stack[i-1].start+stack[i-1].size-1, (arr)+(stack[i-1].start-1)*size, &r);CHKERRQ(ierr);
662*676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
663*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size);CHKERRQ(ierr);
664*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, l, m, r);CHKERRQ(ierr);
665*676f2a66SJacob Faibussowitsch             } else {
666*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size);CHKERRQ(ierr);
667*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, l, m, r);CHKERRQ(ierr);
668*676f2a66SJacob Faibussowitsch             }
669*676f2a66SJacob Faibussowitsch           }
670*676f2a66SJacob Faibussowitsch           /* Update A with merge */
671*676f2a66SJacob Faibussowitsch           stack[i-2].size += stack[i-1].size;
672*676f2a66SJacob Faibussowitsch           /* Push C up the stack */
673*676f2a66SJacob Faibussowitsch           stack[i-1].start = stack[i].start;
674*676f2a66SJacob Faibussowitsch           stack[i-1].size = stack[i].size;
675*676f2a66SJacob Faibussowitsch         } else {
676*676f2a66SJacob Faibussowitsch           /* merge C into B */
677*676f2a66SJacob Faibussowitsch           mergeBC:
678*676f2a66SJacob Faibussowitsch           /* If B[-1] <= C[0] then... you know the drill */
679*676f2a66SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i].start-1)*size, arr+(stack[i].start)*size) > 0) {
680*676f2a66SJacob Faibussowitsch             m = stack[i].start;
681*676f2a66SJacob Faibussowitsch             /* Search B for C[0] insertion */
682*676f2a66SJacob Faibussowitsch             ierr = PetscGallopSearchLeft_Private(arr, size, cmp, stack[i-1].start, stack[i].start-1, arr+stack[i].start*size, &l);CHKERRQ(ierr);
683*676f2a66SJacob Faibussowitsch             /* Search C for B[-1] insertion */
684*676f2a66SJacob Faibussowitsch             ierr = PetscGallopSearchRight_Private(arr, size, cmp, stack[i].start, stack[i].start+stack[i].size-1, (arr)+(stack[i].start-1)*size, &r);CHKERRQ(ierr);
685*676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
686*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(buff, (m-l+1)*size);CHKERRQ(ierr);
687*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortMergeLo_Private(arr, buff->ptr, size, cmp, l, m, r);CHKERRQ(ierr);
688*676f2a66SJacob Faibussowitsch             } else {
689*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(buff, (r-m+1)*size);CHKERRQ(ierr);
690*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortMergeHi_Private(arr, buff->ptr, size, cmp, l, m, r);CHKERRQ(ierr);
691*676f2a66SJacob Faibussowitsch             }
692*676f2a66SJacob Faibussowitsch           }
693*676f2a66SJacob Faibussowitsch           /* Update B with merge */
694*676f2a66SJacob Faibussowitsch           stack[i-1].size += stack[i].size;
695*676f2a66SJacob Faibussowitsch         }
696*676f2a66SJacob Faibussowitsch         --i;
697*676f2a66SJacob Faibussowitsch       } else if (stack[i-1].size <= stack[i].size) {
698*676f2a66SJacob Faibussowitsch         /* merge C into B */
699*676f2a66SJacob Faibussowitsch         goto mergeBC;
700*676f2a66SJacob Faibussowitsch       }
701*676f2a66SJacob Faibussowitsch     }
702*676f2a66SJacob Faibussowitsch     if (itemp == i) break;
703*676f2a66SJacob Faibussowitsch   }
704*676f2a66SJacob Faibussowitsch   *stacksize = i;
705*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
706*676f2a66SJacob Faibussowitsch }
707*676f2a66SJacob Faibussowitsch 
708*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortMergeCollapseWithArray_Private(char *arr, size_t asize, char *barr, size_t bsize, CompFunc cmp, PetscTimSortBuffer *abuff, PetscTimSortBuffer *bbuff, PetscTimSortStack *stack, PetscInt *stacksize)
709*676f2a66SJacob Faibussowitsch {
710*676f2a66SJacob Faibussowitsch   PetscInt       i = *stacksize;
711*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
712*676f2a66SJacob Faibussowitsch 
713*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
714*676f2a66SJacob Faibussowitsch   while (i) {
715*676f2a66SJacob Faibussowitsch     PetscInt l, m, r, itemp = i;
716*676f2a66SJacob Faibussowitsch 
717*676f2a66SJacob Faibussowitsch     if (i == 1) {
718*676f2a66SJacob Faibussowitsch       /* A = stack[i-1], B = stack[i] */
719*676f2a66SJacob Faibussowitsch       if (stack[i-1].size < stack[i].size) {
720*676f2a66SJacob Faibussowitsch         /* if A[-1] <= B[0] then sorted */
721*676f2a66SJacob Faibussowitsch         if ((*cmp)(arr+(stack[i].start-1)*asize, arr+(stack[i].start)*asize) > 0) {
722*676f2a66SJacob Faibussowitsch           m = stack[i].start;
723*676f2a66SJacob Faibussowitsch           /* Search A for B[0] insertion */
724*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, stack[i-1].start, stack[i].start-1, (arr)+stack[i].start*asize, &l);CHKERRQ(ierr);
725*676f2a66SJacob Faibussowitsch           /* Search B for A[-1] insertion */
726*676f2a66SJacob Faibussowitsch           ierr = PetscGallopSearchRight_Private(arr, asize, cmp, stack[i].start, stack[i].start+stack[i].size-1, arr+(stack[i].start-1)*asize, &r);CHKERRQ(ierr);
727*676f2a66SJacob Faibussowitsch           if (m-l <= r-m) {
728*676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize);CHKERRQ(ierr);
729*676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize);CHKERRQ(ierr);
730*676f2a66SJacob Faibussowitsch             ierr = PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, l, m, r);CHKERRQ(ierr);
731*676f2a66SJacob Faibussowitsch           } else {
732*676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize);CHKERRQ(ierr);
733*676f2a66SJacob Faibussowitsch             ierr = PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize);CHKERRQ(ierr);
734*676f2a66SJacob Faibussowitsch             ierr = PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, l, m, r);CHKERRQ(ierr);
735*676f2a66SJacob Faibussowitsch           }
736*676f2a66SJacob Faibussowitsch         }
737*676f2a66SJacob Faibussowitsch         /* Update A with merge */
738*676f2a66SJacob Faibussowitsch         stack[i-1].size += stack[i].size;
739*676f2a66SJacob Faibussowitsch         --i;
740*676f2a66SJacob Faibussowitsch       }
741*676f2a66SJacob Faibussowitsch     } else {
742*676f2a66SJacob Faibussowitsch       /* i > 2, i.e. C exists
743*676f2a66SJacob Faibussowitsch        A = stack[i-2], B = stack[i-1], C = stack[i]; */
744*676f2a66SJacob Faibussowitsch       if (stack[i-2].size <= stack[i-1].size+stack[i].size) {
745*676f2a66SJacob Faibussowitsch         if (stack[i-2].size < stack[i].size) {
746*676f2a66SJacob Faibussowitsch           /* merge B into A, but if A[-1] <= B[0] then already sorted */
747*676f2a66SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i-1].start-1)*asize, arr+(stack[i-1].start)*asize) > 0) {
748*676f2a66SJacob Faibussowitsch             m = stack[i-1].start;
749*676f2a66SJacob Faibussowitsch             /* Search A for B[0] insertion */
750*676f2a66SJacob Faibussowitsch             ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, stack[i-2].start, stack[i-1].start-1, (arr)+(stack[i-1].start)*asize, &l);CHKERRQ(ierr);
751*676f2a66SJacob Faibussowitsch             /* Search B for A[-1] insertion */
752*676f2a66SJacob Faibussowitsch             ierr = PetscGallopSearchRight_Private(arr, asize, cmp, stack[i-1].start, stack[i-1].start+stack[i-1].size-1, (arr)+(stack[i-1].start-1)*asize, &r);CHKERRQ(ierr);
753*676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
754*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize);CHKERRQ(ierr);
755*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize);CHKERRQ(ierr);
756*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, l, m, r);CHKERRQ(ierr);
757*676f2a66SJacob Faibussowitsch             } else {
758*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize);CHKERRQ(ierr);
759*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize);CHKERRQ(ierr);
760*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, l, m, r);CHKERRQ(ierr);
761*676f2a66SJacob Faibussowitsch             }
762*676f2a66SJacob Faibussowitsch           }
763*676f2a66SJacob Faibussowitsch           /* Update A with merge */
764*676f2a66SJacob Faibussowitsch           stack[i-2].size += stack[i-1].size;
765*676f2a66SJacob Faibussowitsch           /* Push C up the stack */
766*676f2a66SJacob Faibussowitsch           stack[i-1].start = stack[i].start;
767*676f2a66SJacob Faibussowitsch           stack[i-1].size = stack[i].size;
768*676f2a66SJacob Faibussowitsch         } else {
769*676f2a66SJacob Faibussowitsch           /* merge C into B */
770*676f2a66SJacob Faibussowitsch           mergeBC:
771*676f2a66SJacob Faibussowitsch           /* If B[-1] <= C[0] then... you know the drill */
772*676f2a66SJacob Faibussowitsch           if ((*cmp)(arr+(stack[i].start-1)*asize, arr+(stack[i].start)*asize) > 0) {
773*676f2a66SJacob Faibussowitsch             m = stack[i].start;
774*676f2a66SJacob Faibussowitsch             /* Search B for C[0] insertion */
775*676f2a66SJacob Faibussowitsch             ierr = PetscGallopSearchLeft_Private(arr, asize, cmp, stack[i-1].start, stack[i].start-1, arr+stack[i].start*asize, &l);CHKERRQ(ierr);
776*676f2a66SJacob Faibussowitsch             /* Search C for B[-1] insertion */
777*676f2a66SJacob Faibussowitsch             ierr = PetscGallopSearchRight_Private(arr, asize, cmp, stack[i].start, stack[i].start+stack[i].size-1, (arr)+(stack[i].start-1)*asize, &r);CHKERRQ(ierr);
778*676f2a66SJacob Faibussowitsch             if (m-l <= r-m) {
779*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(abuff, (m-l+1)*asize);CHKERRQ(ierr);
780*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(bbuff, (m-l+1)*bsize);CHKERRQ(ierr);
781*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortMergeLoWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, l, m, r);CHKERRQ(ierr);
782*676f2a66SJacob Faibussowitsch             } else {
783*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(abuff, (r-m+1)*asize);CHKERRQ(ierr);
784*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortResizeBuffer_Private(bbuff, (r-m+1)*bsize);CHKERRQ(ierr);
785*676f2a66SJacob Faibussowitsch               ierr = PetscTimSortMergeHiWithArray_Private(arr, abuff->ptr, asize, barr, bbuff->ptr, bsize, cmp, l, m, r);CHKERRQ(ierr);
786*676f2a66SJacob Faibussowitsch             }
787*676f2a66SJacob Faibussowitsch           }
788*676f2a66SJacob Faibussowitsch           /* Update B with merge */
789*676f2a66SJacob Faibussowitsch           stack[i-1].size += stack[i].size;
790*676f2a66SJacob Faibussowitsch         }
791*676f2a66SJacob Faibussowitsch         --i;
792*676f2a66SJacob Faibussowitsch       } else if (stack[i-1].size <= stack[i].size) {
793*676f2a66SJacob Faibussowitsch         /* merge C into B */
794*676f2a66SJacob Faibussowitsch         goto mergeBC;
795*676f2a66SJacob Faibussowitsch       }
796*676f2a66SJacob Faibussowitsch     }
797*676f2a66SJacob Faibussowitsch     if (itemp == i) break;
798*676f2a66SJacob Faibussowitsch   }
799*676f2a66SJacob Faibussowitsch   *stacksize = i;
800*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
801*676f2a66SJacob Faibussowitsch }
802*676f2a66SJacob Faibussowitsch 
803*676f2a66SJacob Faibussowitsch /* March sequentially through the array building up a "run" of weakly increasing or strictly decreasing contiguous
804*676f2a66SJacob Faibussowitsch  elements. Decreasing runs are reversed by swapping. If the run is less than minrun, artificially extend it via either
805*676f2a66SJacob Faibussowitsch  binary insertion sort or regulat insertion sort */
806*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortBuildRun_Private(char *arr, char *tarr, size_t size, CompFunc cmp, PetscInt n, PetscInt minrun, PetscInt runstart, PetscInt *runend)
807*676f2a66SJacob Faibussowitsch {
808*676f2a66SJacob Faibussowitsch   const PetscInt re = PetscMin(runstart+minrun, n-1);
809*676f2a66SJacob Faibussowitsch   PetscInt       ri = runstart;
810*676f2a66SJacob Faibussowitsch 
811*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
812*676f2a66SJacob Faibussowitsch   if (PetscUnlikely(runstart == n-1)) {*runend = runstart; PetscFunctionReturn(0);}
813*676f2a66SJacob Faibussowitsch   /* guess whether run is ascending or descending and tally up the longest consecutive run. essentially a coinflip for random data */
814*676f2a66SJacob Faibussowitsch   if ((*cmp)((arr)+(ri+1)*size, (arr)+ri*size) < 0) {
815*676f2a66SJacob Faibussowitsch     ++ri;
816*676f2a66SJacob Faibussowitsch     while (ri < n-1) {
817*676f2a66SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*size, (arr)+ri*size) >= 0) break;
818*676f2a66SJacob Faibussowitsch       ++ri;
819*676f2a66SJacob Faibussowitsch     }
820*676f2a66SJacob Faibussowitsch     {
821*676f2a66SJacob Faibussowitsch       PetscInt lo = runstart, hi = ri;
822*676f2a66SJacob Faibussowitsch       for (; lo < hi; ++lo, --hi) {
823*676f2a66SJacob Faibussowitsch         COPYSWAPPY(arr+lo*size, arr+hi*size, tarr, size);
824*676f2a66SJacob Faibussowitsch       }
825*676f2a66SJacob Faibussowitsch     }
826*676f2a66SJacob Faibussowitsch   } else {
827*676f2a66SJacob Faibussowitsch     ++ri;
828*676f2a66SJacob Faibussowitsch     while (ri < n-1) {
829*676f2a66SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*size, (arr)+ri*size) < 0) break;
830*676f2a66SJacob Faibussowitsch       ++ri;
831*676f2a66SJacob Faibussowitsch     }
832*676f2a66SJacob Faibussowitsch   }
833*676f2a66SJacob Faibussowitsch #if defined(PETSC_USE_DEBUG)
834*676f2a66SJacob Faibussowitsch   {
835*676f2a66SJacob Faibussowitsch     PetscErrorCode ierr;
836*676f2a66SJacob Faibussowitsch     ierr = PetscInfo1(NULL, "natural run length = %D\n", ri-runstart+1);CHKERRQ(ierr);
837*676f2a66SJacob Faibussowitsch   }
838*676f2a66SJacob Faibussowitsch #endif
839*676f2a66SJacob Faibussowitsch   if (ri < re) {
840*676f2a66SJacob Faibussowitsch     /* the attempt failed, this section likely contains random data. If ri got close to minrun (within 50%) then we try
841*676f2a66SJacob Faibussowitsch      binary search */
842*676f2a66SJacob Faibussowitsch     if (ri-runstart <= minrun >> 1) {
843*676f2a66SJacob Faibussowitsch       ++MIN_GALLOP_GLOBAL; /* didn't get close hedge our bets against random data */
844*676f2a66SJacob Faibussowitsch       PetscInsertionSort_Private(arr, tarr, size, cmp, runstart, ri, re);
845*676f2a66SJacob Faibussowitsch     } else {
846*676f2a66SJacob Faibussowitsch       PetscBinaryInsertionSort_Private(arr, tarr, size, cmp, runstart, ri, re);
847*676f2a66SJacob Faibussowitsch     }
848*676f2a66SJacob Faibussowitsch     *runend = re;
849*676f2a66SJacob Faibussowitsch   } else *runend = ri;
850*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
851*676f2a66SJacob Faibussowitsch }
852*676f2a66SJacob Faibussowitsch 
853*676f2a66SJacob Faibussowitsch PETSC_STATIC_INLINE PetscErrorCode PetscTimSortBuildRunWithArray_Private(char *arr, char *atarr, size_t asize, char *barr, char *btarr, size_t bsize, CompFunc cmp, PetscInt n, PetscInt minrun, PetscInt runstart, PetscInt *runend)
854*676f2a66SJacob Faibussowitsch {
855*676f2a66SJacob Faibussowitsch   const PetscInt re = PetscMin(runstart+minrun, n-1);
856*676f2a66SJacob Faibussowitsch   PetscInt       ri = runstart;
857*676f2a66SJacob Faibussowitsch 
858*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
859*676f2a66SJacob Faibussowitsch   if (PetscUnlikely(runstart == n-1)) {*runend = runstart; PetscFunctionReturn(0);}
860*676f2a66SJacob Faibussowitsch   /* guess whether run is ascending or descending and tally up the longest consecutive run. essentially a coinflip for random data */
861*676f2a66SJacob Faibussowitsch   if ((*cmp)((arr)+(ri+1)*asize, arr+(ri*asize)) < 0) {
862*676f2a66SJacob Faibussowitsch     ++ri;
863*676f2a66SJacob Faibussowitsch     while (ri < n-1) {
864*676f2a66SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*asize, (arr)+ri*asize) >= 0) break;
865*676f2a66SJacob Faibussowitsch       ++ri;
866*676f2a66SJacob Faibussowitsch     }
867*676f2a66SJacob Faibussowitsch     {
868*676f2a66SJacob Faibussowitsch       PetscInt lo = runstart, hi = ri;
869*676f2a66SJacob Faibussowitsch       for (; lo < hi; ++lo, --hi) {
870*676f2a66SJacob Faibussowitsch         COPYSWAPPY2(arr+lo*asize, arr+hi*asize, asize, barr+lo*bsize, barr+hi*bsize, bsize, atarr);
871*676f2a66SJacob Faibussowitsch       }
872*676f2a66SJacob Faibussowitsch     }
873*676f2a66SJacob Faibussowitsch   } else {
874*676f2a66SJacob Faibussowitsch     ++ri;
875*676f2a66SJacob Faibussowitsch     while (ri < n-1) {
876*676f2a66SJacob Faibussowitsch       if ((*cmp)((arr)+(ri+1)*asize, (arr)+ri*asize) < 0) break;
877*676f2a66SJacob Faibussowitsch       ++ri;
878*676f2a66SJacob Faibussowitsch     }
879*676f2a66SJacob Faibussowitsch   }
880*676f2a66SJacob Faibussowitsch #if defined(PETSC_USE_DEBUG)
881*676f2a66SJacob Faibussowitsch   {
882*676f2a66SJacob Faibussowitsch     PetscErrorCode ierr;
883*676f2a66SJacob Faibussowitsch     ierr = PetscInfo1(NULL, "natural run length = %D\n", ri-runstart+1);CHKERRQ(ierr);
884*676f2a66SJacob Faibussowitsch   }
885*676f2a66SJacob Faibussowitsch #endif
886*676f2a66SJacob Faibussowitsch   if (ri < re) {
887*676f2a66SJacob Faibussowitsch     /* the attempt failed, this section likely contains random data. If ri got close to minrun (within 50%) then we try
888*676f2a66SJacob Faibussowitsch      binary search */
889*676f2a66SJacob Faibussowitsch     if (ri-runstart <= minrun >> 1) {
890*676f2a66SJacob Faibussowitsch       ++MIN_GALLOP_GLOBAL; /* didn't get close hedge our bets against random data */
891*676f2a66SJacob Faibussowitsch       PetscInsertionSortWithArray_Private(arr, atarr, asize, barr, btarr, bsize, cmp, runstart, ri, re);
892*676f2a66SJacob Faibussowitsch     } else {
893*676f2a66SJacob Faibussowitsch       PetscBinaryInsertionSortWithArray_Private(arr, atarr, asize, barr, btarr, bsize, cmp, runstart, ri, re);
894*676f2a66SJacob Faibussowitsch     }
895*676f2a66SJacob Faibussowitsch     *runend = re;
896*676f2a66SJacob Faibussowitsch   } else *runend = ri;
897*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
898*676f2a66SJacob Faibussowitsch }
899*676f2a66SJacob Faibussowitsch 
900*676f2a66SJacob Faibussowitsch /*
901*676f2a66SJacob Faibussowitsch   PetscTimSort - Sorts an array in place in increasing order using Tim Peters adaptive sorting algorithm.
902*676f2a66SJacob Faibussowitsch 
903*676f2a66SJacob Faibussowitsch   Not Collective
904*676f2a66SJacob Faibussowitsch 
905*676f2a66SJacob Faibussowitsch   Input Parameters:
906*676f2a66SJacob Faibussowitsch + n    - number of values
907*676f2a66SJacob Faibussowitsch . arr  - array to be sorted
908*676f2a66SJacob Faibussowitsch . size - size in bytes of the datatype held in arr
909*676f2a66SJacob Faibussowitsch - cmp  - function pointer to comparison function
910*676f2a66SJacob Faibussowitsch 
911*676f2a66SJacob Faibussowitsch   Output Parameters:
912*676f2a66SJacob Faibussowitsch . arr  - sorted array
913*676f2a66SJacob Faibussowitsch 
914*676f2a66SJacob Faibussowitsch   Sample usage:
915*676f2a66SJacob Faibussowitsch   The comparison function must follow the qsort() comparison function paradigm, returning the sign of the difference
916*676f2a66SJacob Faibussowitsch   between its arguments. If left < right : return -1, if left == right : return 0, if left > right : return 1. The user
917*676f2a66SJacob Faibussowitsch   may also
918*676f2a66SJacob Faibussowitsch  change or reverse the order of the sort by flipping the above. Note that stability of the sort is only guaranteed if
919*676f2a66SJacob Faibussowitsch  the comparison function forms a valid trigraph. For example when sorting an array of type "my_type" in increasing
920*676f2a66SJacob Faibussowitsch   order:
921*676f2a66SJacob Faibussowitsch 
922*676f2a66SJacob Faibussowitsch .vb
923*676f2a66SJacob Faibussowitsch   int my_increasing_comparison_function(const void *left, const void *right) {
924*676f2a66SJacob Faibussowitsch     my_type l = *(my_type *) left, r = *(my_type *) right;
925*676f2a66SJacob Faibussowitsch     return (l < r) ? -1 : (l > r);
926*676f2a66SJacob Faibussowitsch   }
927*676f2a66SJacob Faibussowitsch .ve
928*676f2a66SJacob Faibussowitsch   Then pass the function
929*676f2a66SJacob Faibussowitsch .vb
930*676f2a66SJacob Faibussowitsch   PetscTimSort(n, arr, sizeof(arr[0]), my_increasing_comparison_function)
931*676f2a66SJacob Faibussowitsch .ve
932*676f2a66SJacob Faibussowitsch 
933*676f2a66SJacob Faibussowitsch   Notes: Timsort makes the assumption that input data is already likely partially ordered, or that it contains
934*676f2a66SJacob Faibussowitsch   contiguous sections (termed 'runs') where the data is locally ordered (but not necessarily globally ordered). It
935*676f2a66SJacob Faibussowitsch   therefore aims
936*676f2a66SJacob Faibussowitsch  to select slices of the array in such a way that resulting mergesorts operate on near perfectly length-balanced
937*676f2a66SJacob Faibussowitsch  arrays. To do so it repeatedly triggers attempts throughout to merge adjacent runs.
938*676f2a66SJacob Faibussowitsch 
939*676f2a66SJacob Faibussowitsch   Should one run continuously "win" a comparison the algorithm begins the "gallop" phase. It will aggressively search
940*676f2a66SJacob Faibussowitsch   the "winner" for the location of the "losers" next entry (and vice versa) to copy all preceding elements into place in
941*676f2a66SJacob Faibussowitsch   bulk. However if the data is truly unordered (as is the case with random data) the immense gains possible from these
942*676f2a66SJacob Faibussowitsch   searches are expected __not__ to repay their costs. While adjacent arrays are almost all nearly the same size, they
943*676f2a66SJacob Faibussowitsch   likely all contain similar data.
944*676f2a66SJacob Faibussowitsch 
945*676f2a66SJacob Faibussowitsch   References:
946*676f2a66SJacob Faibussowitsch   1. - Tim Peters. https://bugs.python.org/file4451/timsort.txt
947*676f2a66SJacob Faibussowitsch 
948*676f2a66SJacob Faibussowitsch   Level: developer
949*676f2a66SJacob Faibussowitsch 
950*676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscIntSortSemiOrdered(), PetscRealSortSemiOrdered(), PetscMPIIntSortSemiOrdered()
951*676f2a66SJacob Faibussowitsch */
952*676f2a66SJacob Faibussowitsch PetscErrorCode PetscTimSort(PetscInt n, void *arr, size_t size, int (*cmp)(const void *, const void *))
953*676f2a66SJacob Faibussowitsch {
954*676f2a66SJacob Faibussowitsch   PetscInt           stacksize = 0, minrun, runstart = 0, runend = 0;
955*676f2a66SJacob Faibussowitsch   PetscTimSortStack  runstack[128];
956*676f2a66SJacob Faibussowitsch   PetscTimSortBuffer buff;
957*676f2a66SJacob Faibussowitsch   PetscErrorCode     ierr;
958*676f2a66SJacob Faibussowitsch   /* stacksize  = log_phi(n) = log_2(n)/log_2(phi), so 128 is enough for ~5.614e26 elements.
959*676f2a66SJacob Faibussowitsch    It is so unlikely that this limit is reached that this is __never__ checked for */
960*676f2a66SJacob Faibussowitsch 
961*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
962*676f2a66SJacob Faibussowitsch   /* Compute minrun. Minrun should be (32, 65) such that N/minrun
963*676f2a66SJacob Faibussowitsch    is a power of 2 or one plus a power of 2 */
964*676f2a66SJacob Faibussowitsch   {
965*676f2a66SJacob Faibussowitsch     PetscInt t = n, r = 0;
966*676f2a66SJacob Faibussowitsch     /* r becomes 1 if the least significant bits contain at least one off bit */
967*676f2a66SJacob Faibussowitsch     while (t >= 64) {
968*676f2a66SJacob Faibussowitsch       r |= t & 1;
969*676f2a66SJacob Faibussowitsch       t >>= 1;
970*676f2a66SJacob Faibussowitsch     }
971*676f2a66SJacob Faibussowitsch     minrun = t + r;
972*676f2a66SJacob Faibussowitsch   }
973*676f2a66SJacob Faibussowitsch   if (PetscUnlikelyDebug(minrun < 32 || minrun > 65)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Calculated minrun %D not in range (32,65)",minrun);
974*676f2a66SJacob Faibussowitsch   ierr = PetscInfo1(NULL, "minrun = %D\n", minrun);CHKERRQ(ierr);
975*676f2a66SJacob Faibussowitsch   ierr = PetscMalloc1((size_t) minrun*size, &buff.ptr);CHKERRQ(ierr);
976*676f2a66SJacob Faibussowitsch   buff.size = (size_t) minrun*size;
977*676f2a66SJacob Faibussowitsch   buff.maxsize = (size_t) n*size;
978*676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
979*676f2a66SJacob Faibussowitsch   while (runstart < n) {
980*676f2a66SJacob Faibussowitsch     /* Check if additional entries are at least partially ordered and build natural run */
981*676f2a66SJacob Faibussowitsch     ierr = PetscTimSortBuildRun_Private((char *)arr, buff.ptr, size, cmp, n, minrun, runstart, &runend);CHKERRQ(ierr);
982*676f2a66SJacob Faibussowitsch     runstack[stacksize].start = runstart;
983*676f2a66SJacob Faibussowitsch     runstack[stacksize].size = runend-runstart+1;
984*676f2a66SJacob Faibussowitsch     ierr = PetscTimSortMergeCollapse_Private((char *)arr, size, cmp, &buff, runstack, &stacksize);CHKERRQ(ierr);
985*676f2a66SJacob Faibussowitsch     ++stacksize;
986*676f2a66SJacob Faibussowitsch     runstart = runend+1;
987*676f2a66SJacob Faibussowitsch   }
988*676f2a66SJacob Faibussowitsch   /* Have been inside while, so discard last stacksize++ */
989*676f2a66SJacob Faibussowitsch   --stacksize;
990*676f2a66SJacob Faibussowitsch   ierr = PetscTimSortForceCollapse_Private((char *)arr, size, cmp, &buff, runstack, stacksize);CHKERRQ(ierr);
991*676f2a66SJacob Faibussowitsch   ierr = PetscFree(buff.ptr);CHKERRQ(ierr);
992*676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
993*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
994*676f2a66SJacob Faibussowitsch }
995*676f2a66SJacob Faibussowitsch 
996*676f2a66SJacob Faibussowitsch /*
997*676f2a66SJacob Faibussowitsch   PetscTimSortWithArray - Sorts an array in place in increasing order using Tim Peters adaptive sorting algorithm and
998*676f2a66SJacob Faibussowitsch   reorders a second array to match the first. The arrays need not be the same type.
999*676f2a66SJacob Faibussowitsch 
1000*676f2a66SJacob Faibussowitsch   Not Collective
1001*676f2a66SJacob Faibussowitsch 
1002*676f2a66SJacob Faibussowitsch   Input Parameters:
1003*676f2a66SJacob Faibussowitsch + n     - number of values
1004*676f2a66SJacob Faibussowitsch . arr   - array to be sorted
1005*676f2a66SJacob Faibussowitsch . asize - size in bytes of the datatype held in arr
1006*676f2a66SJacob Faibussowitsch . barr  - array to be reordered
1007*676f2a66SJacob Faibussowitsch . asize - size in bytes of the datatype held in barr
1008*676f2a66SJacob Faibussowitsch - cmp   - function pointer to comparison function
1009*676f2a66SJacob Faibussowitsch 
1010*676f2a66SJacob Faibussowitsch   Output Parameters:
1011*676f2a66SJacob Faibussowitsch + arr  - sorted array
1012*676f2a66SJacob Faibussowitsch + barr - reordered array
1013*676f2a66SJacob Faibussowitsch 
1014*676f2a66SJacob Faibussowitsch   Sample usage:
1015*676f2a66SJacob Faibussowitsch   The comparison function must follow the qsort() comparison function paradigm, returning the sign of the difference
1016*676f2a66SJacob Faibussowitsch   between its arguments. If left < right : return -1, if left == right : return 0, if left > right : return 1. The user
1017*676f2a66SJacob Faibussowitsch   may also change or reverse the order of the sort by flipping the above. Note that stability of the sort is only
1018*676f2a66SJacob Faibussowitsch   guaranteed if the comparison function forms a valid trigraph. For example when sorting an array of type "my_type" in
1019*676f2a66SJacob Faibussowitsch   increasing order:
1020*676f2a66SJacob Faibussowitsch 
1021*676f2a66SJacob Faibussowitsch .vb
1022*676f2a66SJacob Faibussowitsch   int my_increasing_comparison_function(const void *left, const void *right) {
1023*676f2a66SJacob Faibussowitsch     my_type l = *(my_type *) left, r = *(my_type *) right;
1024*676f2a66SJacob Faibussowitsch     return (l < r) ? -1 : (l > r);
1025*676f2a66SJacob Faibussowitsch   }
1026*676f2a66SJacob Faibussowitsch .ve
1027*676f2a66SJacob Faibussowitsch   Then pass the function
1028*676f2a66SJacob Faibussowitsch .vb
1029*676f2a66SJacob Faibussowitsch   PetscTimSortWithArray(n, arr, sizeof(arr[0]), barr, sizeof(barr[0]), my_increasing_comparison_function)
1030*676f2a66SJacob Faibussowitsch .ve
1031*676f2a66SJacob Faibussowitsch 
1032*676f2a66SJacob Faibussowitsch   Notes:
1033*676f2a66SJacob 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
1034*676f2a66SJacob Faibussowitsch   overlap.
1035*676f2a66SJacob Faibussowitsch 
1036*676f2a66SJacob Faibussowitsch   Timsort makes the assumption that input data is already likely partially ordered, or that it contains contiguous
1037*676f2a66SJacob Faibussowitsch   sections (termed 'runs') where the data is locally ordered (but not necessarily globally ordered). It therefore aims
1038*676f2a66SJacob Faibussowitsch  to select slices of the array in such a way that resulting mergesorts operate on near perfectly length-balanced
1039*676f2a66SJacob Faibussowitsch  arrays. To do so it repeatedly triggers attempts throughout to merge adjacent runs.
1040*676f2a66SJacob Faibussowitsch 
1041*676f2a66SJacob Faibussowitsch   Should one run continuously "win" a comparison the algorithm begins the "gallop" phase. It will aggressively
1042*676f2a66SJacob Faibussowitsch   search the "winner" for the location of the "losers" next entry (and vice versa) to copy all preceding elements into
1043*676f2a66SJacob Faibussowitsch   place in bulk. However if the data is truly unordered (as is the case with random data) the immense gains possible
1044*676f2a66SJacob Faibussowitsch   from these searches are expected __not__ to repay their costs. While adjacent arrays are almost all nearly the same
1045*676f2a66SJacob Faibussowitsch   size, they likely all contain similar data.
1046*676f2a66SJacob Faibussowitsch 
1047*676f2a66SJacob Faibussowitsch   References:
1048*676f2a66SJacob Faibussowitsch   1. - Tim Peters. https://bugs.python.org/file4451/timsort.txt
1049*676f2a66SJacob Faibussowitsch 
1050*676f2a66SJacob Faibussowitsch   Level: developer
1051*676f2a66SJacob Faibussowitsch 
1052*676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscIntSortSemiOrderedWithArray(), PetscRealSortSemiOrderedWithArrayInt(), PetscMPIIntSortSemiOrderedWithArray()
1053*676f2a66SJacob Faibussowitsch */
1054*676f2a66SJacob Faibussowitsch PetscErrorCode PetscTimSortWithArray(PetscInt n, void *arr, size_t asize, void *barr, size_t bsize, int (*cmp)(const void *, const void *))
1055*676f2a66SJacob Faibussowitsch {
1056*676f2a66SJacob Faibussowitsch   PetscInt           stacksize = 0, minrun, runstart = 0, runend = 0;
1057*676f2a66SJacob Faibussowitsch   PetscTimSortStack  runstack[128];
1058*676f2a66SJacob Faibussowitsch   PetscTimSortBuffer abuff, bbuff;
1059*676f2a66SJacob Faibussowitsch   PetscErrorCode     ierr;
1060*676f2a66SJacob Faibussowitsch   /* stacksize  = log_phi(n) = log_2(n)/log_2(phi), so 128 is enough for ~5.614e26 elements.
1061*676f2a66SJacob Faibussowitsch    It is so unlikely that this limit is reached that this is __never__ checked for */
1062*676f2a66SJacob Faibussowitsch 
1063*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1064*676f2a66SJacob Faibussowitsch   /* Compute minrun. Minrun should be (32, 65) such that N/minrun
1065*676f2a66SJacob Faibussowitsch    is a power of 2 or one plus a power of 2 */
1066*676f2a66SJacob Faibussowitsch   {
1067*676f2a66SJacob Faibussowitsch     PetscInt t = n, r = 0;
1068*676f2a66SJacob Faibussowitsch     /* r becomes 1 if the least significant bits contain at least one off bit */
1069*676f2a66SJacob Faibussowitsch     while (t >= 64) {
1070*676f2a66SJacob Faibussowitsch       r |= t & 1;
1071*676f2a66SJacob Faibussowitsch       t >>= 1;
1072*676f2a66SJacob Faibussowitsch     }
1073*676f2a66SJacob Faibussowitsch     minrun = t + r;
1074*676f2a66SJacob Faibussowitsch   }
1075*676f2a66SJacob Faibussowitsch   if (PetscUnlikelyDebug(minrun < 32 || minrun > 65)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Calculated minrun %D not in range (32,65)",minrun);
1076*676f2a66SJacob Faibussowitsch   ierr = PetscInfo1(NULL, "minrun = %D\n", minrun);CHKERRQ(ierr);
1077*676f2a66SJacob Faibussowitsch   ierr = PetscMalloc1((size_t) minrun*asize, &abuff.ptr);CHKERRQ(ierr);
1078*676f2a66SJacob Faibussowitsch   abuff.size = (size_t) minrun*asize;
1079*676f2a66SJacob Faibussowitsch   abuff.maxsize = (size_t) n*asize;
1080*676f2a66SJacob Faibussowitsch   ierr = PetscMalloc1((size_t) minrun*bsize, &bbuff.ptr);CHKERRQ(ierr);
1081*676f2a66SJacob Faibussowitsch   bbuff.size = (size_t) minrun*bsize;
1082*676f2a66SJacob Faibussowitsch   bbuff.maxsize = (size_t) n*bsize;
1083*676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
1084*676f2a66SJacob Faibussowitsch   while (runstart < n) {
1085*676f2a66SJacob Faibussowitsch     /* Check if additional entries are at least partially ordered and build natural run */
1086*676f2a66SJacob Faibussowitsch     ierr = PetscTimSortBuildRunWithArray_Private((char *)arr, abuff.ptr, asize, (char *)barr, bbuff.ptr, bsize, cmp, n, minrun, runstart, &runend);CHKERRQ(ierr);
1087*676f2a66SJacob Faibussowitsch     runstack[stacksize].start = runstart;
1088*676f2a66SJacob Faibussowitsch     runstack[stacksize].size = runend-runstart+1;
1089*676f2a66SJacob Faibussowitsch     ierr = PetscTimSortMergeCollapseWithArray_Private((char *)arr, asize, (char *)barr, bsize, cmp, &abuff, &bbuff, runstack, &stacksize);CHKERRQ(ierr);
1090*676f2a66SJacob Faibussowitsch     ++stacksize;
1091*676f2a66SJacob Faibussowitsch     runstart = runend+1;
1092*676f2a66SJacob Faibussowitsch   }
1093*676f2a66SJacob Faibussowitsch   /* Have been inside while, so discard last stacksize++ */
1094*676f2a66SJacob Faibussowitsch   --stacksize;
1095*676f2a66SJacob Faibussowitsch   ierr = PetscTimSortForceCollapseWithArray_Private((char *)arr, asize, (char *)barr, bsize, cmp, &abuff, &bbuff, runstack, stacksize);CHKERRQ(ierr);
1096*676f2a66SJacob Faibussowitsch   ierr = PetscFree(abuff.ptr);CHKERRQ(ierr);
1097*676f2a66SJacob Faibussowitsch   ierr = PetscFree(bbuff.ptr);CHKERRQ(ierr);
1098*676f2a66SJacob Faibussowitsch   MIN_GALLOP_GLOBAL = MIN_GALLOP_CONST_GLOBAL;
1099*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1100*676f2a66SJacob Faibussowitsch }
1101*676f2a66SJacob Faibussowitsch 
1102*676f2a66SJacob Faibussowitsch /*@
1103*676f2a66SJacob Faibussowitsch    PetscIntSortSemiOrdered - Sorts an array of integers in place in increasing order.
1104*676f2a66SJacob Faibussowitsch 
1105*676f2a66SJacob Faibussowitsch    Not Collective
1106*676f2a66SJacob Faibussowitsch 
1107*676f2a66SJacob Faibussowitsch    Input Parameters:
1108*676f2a66SJacob Faibussowitsch +  n   - number of values
1109*676f2a66SJacob Faibussowitsch -  arr - array of integers
1110*676f2a66SJacob Faibussowitsch 
1111*676f2a66SJacob Faibussowitsch    Output Parameters:
1112*676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1113*676f2a66SJacob Faibussowitsch 
1114*676f2a66SJacob Faibussowitsch    Notes:
1115*676f2a66SJacob Faibussowitsch    If the array is less than 64 entries long PetscSortInt() is automatically used.
1116*676f2a66SJacob Faibussowitsch 
1117*676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortInt(). While this function works for any array of integers it is
1118*676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1119*676f2a66SJacob Faibussowitsch    recomended that the user benchmark their code to see which routine is fastest.
1120*676f2a66SJacob Faibussowitsch 
1121*676f2a66SJacob Faibussowitsch    Level: intermediate
1122*676f2a66SJacob Faibussowitsch 
1123*676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscSortInt(), PetscSortIntWithPermutation()
1124*676f2a66SJacob Faibussowitsch @*/
1125*676f2a66SJacob Faibussowitsch PetscErrorCode PetscIntSortSemiOrdered(PetscInt n, PetscInt arr[])
1126*676f2a66SJacob Faibussowitsch {
1127*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
1128*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1129*676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1130*676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr,2);
1131*676f2a66SJacob Faibussowitsch   if (n < 64) {
1132*676f2a66SJacob Faibussowitsch     ierr = PetscSortInt(n, arr);CHKERRQ(ierr);
1133*676f2a66SJacob Faibussowitsch   } else {
1134*676f2a66SJacob Faibussowitsch     ierr = PetscTimSort(n, arr, sizeof(PetscInt), Compare_PetscInt_Private);CHKERRQ(ierr);
1135*676f2a66SJacob Faibussowitsch   }
1136*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1137*676f2a66SJacob Faibussowitsch }
1138*676f2a66SJacob Faibussowitsch 
1139*676f2a66SJacob Faibussowitsch /*@
1140*676f2a66SJacob Faibussowitsch    PetscIntSortSemiOrderedWithArray - Sorts an array of integers in place in increasing order and reorders a second
1141*676f2a66SJacob Faibussowitsch    array to match the first.
1142*676f2a66SJacob Faibussowitsch 
1143*676f2a66SJacob Faibussowitsch    Not Collective
1144*676f2a66SJacob Faibussowitsch 
1145*676f2a66SJacob Faibussowitsch    Input Parameters:
1146*676f2a66SJacob Faibussowitsch +  n   - number of values
1147*676f2a66SJacob Faibussowitsch .  arr1 - array of integers to be sorted
1148*676f2a66SJacob Faibussowitsch -  arr2 - array of integers to be reordered
1149*676f2a66SJacob Faibussowitsch 
1150*676f2a66SJacob Faibussowitsch    Output Parameters:
1151*676f2a66SJacob Faibussowitsch +  arr1 - sorted array of integers
1152*676f2a66SJacob Faibussowitsch -  arr2 - reordered array of integers
1153*676f2a66SJacob Faibussowitsch 
1154*676f2a66SJacob Faibussowitsch    Notes:
1155*676f2a66SJacob Faibussowitsch    The arrays CANNOT overlap.
1156*676f2a66SJacob Faibussowitsch 
1157*676f2a66SJacob Faibussowitsch    If the array to be sorted is less than 64 entries long PetscSortIntWithArray() is automatically used.
1158*676f2a66SJacob Faibussowitsch 
1159*676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortIntWithArray(). While this function works for any array of integers it is
1160*676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1161*676f2a66SJacob Faibussowitsch    recomended that the user benchmark their code to see which routine is fastest.
1162*676f2a66SJacob Faibussowitsch 
1163*676f2a66SJacob Faibussowitsch    Level: intermediate
1164*676f2a66SJacob Faibussowitsch 
1165*676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscSortIntWithArray(), PetscSortIntWithPermutation()
1166*676f2a66SJacob Faibussowitsch @*/
1167*676f2a66SJacob Faibussowitsch PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt n, PetscInt arr1[], PetscInt arr2[])
1168*676f2a66SJacob Faibussowitsch {
1169*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
1170*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1171*676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr1,2);
1172*676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2,3);
1173*676f2a66SJacob Faibussowitsch   if (n == 1) PetscFunctionReturn(0);
1174*676f2a66SJacob Faibussowitsch   if (n < 64) {
1175*676f2a66SJacob Faibussowitsch     ierr = PetscSortIntWithArray(n, arr1, arr2);CHKERRQ(ierr);
1176*676f2a66SJacob Faibussowitsch   } else {
1177*676f2a66SJacob Faibussowitsch     ierr = PetscTimSortWithArray(n, arr1, sizeof(PetscInt), arr2, sizeof(PetscInt), Compare_PetscInt_Private);CHKERRQ(ierr);
1178*676f2a66SJacob Faibussowitsch   }
1179*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1180*676f2a66SJacob Faibussowitsch }
1181*676f2a66SJacob Faibussowitsch 
1182*676f2a66SJacob Faibussowitsch /*@
1183*676f2a66SJacob Faibussowitsch    PetscMPIIntSortSemiOrdered - Sorts an array of PetscMPIInts in place in increasing order.
1184*676f2a66SJacob Faibussowitsch 
1185*676f2a66SJacob Faibussowitsch    Not Collective
1186*676f2a66SJacob Faibussowitsch 
1187*676f2a66SJacob Faibussowitsch    Input Parameters:
1188*676f2a66SJacob Faibussowitsch +  n   - number of values
1189*676f2a66SJacob Faibussowitsch -  arr - array of PetscMPIInts
1190*676f2a66SJacob Faibussowitsch 
1191*676f2a66SJacob Faibussowitsch    Output Parameters:
1192*676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1193*676f2a66SJacob Faibussowitsch 
1194*676f2a66SJacob Faibussowitsch    Notes:
1195*676f2a66SJacob Faibussowitsch    If the array is less than 64 entries long PetscSortMPIInt() is automatically used.
1196*676f2a66SJacob Faibussowitsch 
1197*676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortMPIInt(). While this function works for any array of PetscMPIInts it is
1198*676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1199*676f2a66SJacob Faibussowitsch    recomended that the user benchmark their code to see which routine is fastest.
1200*676f2a66SJacob Faibussowitsch 
1201*676f2a66SJacob Faibussowitsch    Level: intermediate
1202*676f2a66SJacob Faibussowitsch 
1203*676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscSortMPIInt()
1204*676f2a66SJacob Faibussowitsch @*/
1205*676f2a66SJacob Faibussowitsch PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt n, PetscMPIInt arr[])
1206*676f2a66SJacob Faibussowitsch {
1207*676f2a66SJacob Faibussowitsch   PetscErrorCode  ierr;
1208*676f2a66SJacob Faibussowitsch 
1209*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1210*676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr,2);
1211*676f2a66SJacob Faibussowitsch   if (n == 1) PetscFunctionReturn(0);
1212*676f2a66SJacob Faibussowitsch   if (n < 64) {
1213*676f2a66SJacob Faibussowitsch     ierr = PetscSortMPIInt(n, arr);CHKERRQ(ierr);
1214*676f2a66SJacob Faibussowitsch   } else {
1215*676f2a66SJacob Faibussowitsch     ierr = PetscTimSort(n, arr, sizeof(PetscMPIInt), Compare_PetscMPIInt_Private);CHKERRQ(ierr);
1216*676f2a66SJacob Faibussowitsch   }
1217*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1218*676f2a66SJacob Faibussowitsch }
1219*676f2a66SJacob Faibussowitsch 
1220*676f2a66SJacob Faibussowitsch /*@
1221*676f2a66SJacob Faibussowitsch    PetscMPIIntSortSemiOrderedWithArray - Sorts an array of integers in place in increasing order and reorders a second
1222*676f2a66SJacob Faibussowitsch    array to match the first.
1223*676f2a66SJacob Faibussowitsch 
1224*676f2a66SJacob Faibussowitsch    Not Collective
1225*676f2a66SJacob Faibussowitsch 
1226*676f2a66SJacob Faibussowitsch    Input Parameters:
1227*676f2a66SJacob Faibussowitsch +  n   - number of values
1228*676f2a66SJacob Faibussowitsch .  arr1 - array of integers to be sorted
1229*676f2a66SJacob Faibussowitsch -  arr2 - array of integers to be reordered
1230*676f2a66SJacob Faibussowitsch 
1231*676f2a66SJacob Faibussowitsch    Output Parameters:
1232*676f2a66SJacob Faibussowitsch +  arr1 - sorted array of integers
1233*676f2a66SJacob Faibussowitsch -  arr2 - reordered array of integers
1234*676f2a66SJacob Faibussowitsch 
1235*676f2a66SJacob Faibussowitsch    Notes:
1236*676f2a66SJacob Faibussowitsch    The arrays CANNOT overlap.
1237*676f2a66SJacob Faibussowitsch 
1238*676f2a66SJacob Faibussowitsch    If the array to be sorted is less than 64 entries long PetscSortMPIIntWithArray() is automatically used.
1239*676f2a66SJacob Faibussowitsch 
1240*676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortMPIIntWithArray(). While this function works for any array of integers it is
1241*676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1242*676f2a66SJacob Faibussowitsch    recomended that the user benchmark their code to see which routine is fastest.
1243*676f2a66SJacob Faibussowitsch 
1244*676f2a66SJacob Faibussowitsch    Level: intermediate
1245*676f2a66SJacob Faibussowitsch 
1246*676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscSortMPIIntWithArray(), PetscSortMPIIntWithPermutation()
1247*676f2a66SJacob Faibussowitsch @*/
1248*676f2a66SJacob Faibussowitsch PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt n, PetscMPIInt arr1[], PetscMPIInt arr2[])
1249*676f2a66SJacob Faibussowitsch {
1250*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
1251*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1252*676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1253*676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr1,2);
1254*676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2,3);
1255*676f2a66SJacob Faibussowitsch   if (n < 64) {
1256*676f2a66SJacob Faibussowitsch     ierr = PetscSortMPIIntWithArray(n, arr1, arr2);CHKERRQ(ierr);
1257*676f2a66SJacob Faibussowitsch   } else {
1258*676f2a66SJacob Faibussowitsch     ierr = PetscTimSortWithArray(n, arr1, sizeof(PetscMPIInt), arr2, sizeof(PetscMPIInt), Compare_PetscMPIInt_Private);CHKERRQ(ierr);
1259*676f2a66SJacob Faibussowitsch   }
1260*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1261*676f2a66SJacob Faibussowitsch }
1262*676f2a66SJacob Faibussowitsch 
1263*676f2a66SJacob Faibussowitsch /*@
1264*676f2a66SJacob Faibussowitsch    PetscRealSortSemiOrdered - Sorts an array of PetscReals in place in increasing order.
1265*676f2a66SJacob Faibussowitsch 
1266*676f2a66SJacob Faibussowitsch    Not Collective
1267*676f2a66SJacob Faibussowitsch 
1268*676f2a66SJacob Faibussowitsch    Input Parameters:
1269*676f2a66SJacob Faibussowitsch +  n   - number of values
1270*676f2a66SJacob Faibussowitsch -  arr - array of PetscReals
1271*676f2a66SJacob Faibussowitsch 
1272*676f2a66SJacob Faibussowitsch    Output Parameters:
1273*676f2a66SJacob Faibussowitsch .  arr - sorted array of integers
1274*676f2a66SJacob Faibussowitsch 
1275*676f2a66SJacob Faibussowitsch    Notes:
1276*676f2a66SJacob Faibussowitsch    If the array is less than 64 entries long PetscSortReal() is automatically used.
1277*676f2a66SJacob Faibussowitsch 
1278*676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortReal(). While this function works for any array of PetscReals it is
1279*676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1280*676f2a66SJacob Faibussowitsch    recomended that the user benchmark their code to see which routine is fastest.
1281*676f2a66SJacob Faibussowitsch 
1282*676f2a66SJacob Faibussowitsch    Level: intermediate
1283*676f2a66SJacob Faibussowitsch 
1284*676f2a66SJacob Faibussowitsch .seealso: PetscTimSort(), PetscSortReal(), PetscSortRealWithPermutation()
1285*676f2a66SJacob Faibussowitsch @*/
1286*676f2a66SJacob Faibussowitsch PetscErrorCode PetscRealSortSemiOrdered(PetscInt n, PetscReal arr[])
1287*676f2a66SJacob Faibussowitsch {
1288*676f2a66SJacob Faibussowitsch   PetscErrorCode  ierr;
1289*676f2a66SJacob Faibussowitsch 
1290*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1291*676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1292*676f2a66SJacob Faibussowitsch   PetscValidRealPointer(arr,2);
1293*676f2a66SJacob Faibussowitsch   if (n < 64) {
1294*676f2a66SJacob Faibussowitsch     ierr = PetscSortReal(n, arr);CHKERRQ(ierr);
1295*676f2a66SJacob Faibussowitsch   } else {
1296*676f2a66SJacob Faibussowitsch     ierr = PetscTimSort(n, arr, sizeof(PetscReal), Compare_PetscReal_Private);CHKERRQ(ierr);
1297*676f2a66SJacob Faibussowitsch   }
1298*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1299*676f2a66SJacob Faibussowitsch }
1300*676f2a66SJacob Faibussowitsch 
1301*676f2a66SJacob Faibussowitsch /*@
1302*676f2a66SJacob Faibussowitsch    PetscRealSortSemiOrderedWithArrayInt - Sorts an array of PetscReals in place in increasing order and reorders a second
1303*676f2a66SJacob Faibussowitsch    array of PetscInts to match the first.
1304*676f2a66SJacob Faibussowitsch 
1305*676f2a66SJacob Faibussowitsch    Not Collective
1306*676f2a66SJacob Faibussowitsch 
1307*676f2a66SJacob Faibussowitsch    Input Parameters:
1308*676f2a66SJacob Faibussowitsch +  n   - number of values
1309*676f2a66SJacob Faibussowitsch .  arr1 - array of PetscReals to be sorted
1310*676f2a66SJacob Faibussowitsch -  arr2 - array of PetscReals to be reordered
1311*676f2a66SJacob Faibussowitsch 
1312*676f2a66SJacob Faibussowitsch    Output Parameters:
1313*676f2a66SJacob Faibussowitsch +  arr1 - sorted array of PetscReals
1314*676f2a66SJacob Faibussowitsch -  arr2 - reordered array of PetscInts
1315*676f2a66SJacob Faibussowitsch 
1316*676f2a66SJacob Faibussowitsch    Notes:
1317*676f2a66SJacob Faibussowitsch    If the array to be sorted is less than 64 entries long PetscSortRealWithArrayInt() is automatically used.
1318*676f2a66SJacob Faibussowitsch 
1319*676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscSortRealWithArray(). While this function works for any array of PetscReals it is
1320*676f2a66SJacob Faibussowitsch    significantly faster if the array is not totally random. There are exceptions to this and so it is __highly__
1321*676f2a66SJacob Faibussowitsch    recomended that the user benchmark their code to see which routine is fastest.
1322*676f2a66SJacob Faibussowitsch 
1323*676f2a66SJacob Faibussowitsch    Level: intermediate
1324*676f2a66SJacob Faibussowitsch 
1325*676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscSortRealWithArrayInt(), PetscSortRealWithPermutation()
1326*676f2a66SJacob Faibussowitsch @*/
1327*676f2a66SJacob Faibussowitsch PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt n, PetscReal arr1[], PetscInt arr2[])
1328*676f2a66SJacob Faibussowitsch {
1329*676f2a66SJacob Faibussowitsch   PetscErrorCode ierr;
1330*676f2a66SJacob Faibussowitsch   PetscFunctionBegin;
1331*676f2a66SJacob Faibussowitsch   if (n <= 1) PetscFunctionReturn(0);
1332*676f2a66SJacob Faibussowitsch   PetscValidRealPointer(arr1,2);
1333*676f2a66SJacob Faibussowitsch   PetscValidIntPointer(arr2,3);
1334*676f2a66SJacob Faibussowitsch   if (n < 64) {
1335*676f2a66SJacob Faibussowitsch     ierr = PetscSortRealWithArrayInt(n, arr1, arr2);CHKERRQ(ierr);
1336*676f2a66SJacob Faibussowitsch   } else {
1337*676f2a66SJacob Faibussowitsch     ierr = PetscTimSortWithArray(n, arr1, sizeof(PetscReal), arr2, sizeof(PetscInt), Compare_PetscReal_Private);CHKERRQ(ierr);
1338*676f2a66SJacob Faibussowitsch   }
1339*676f2a66SJacob Faibussowitsch   PetscFunctionReturn(0);
1340*676f2a66SJacob Faibussowitsch }
1341