xref: /petsc/src/sys/utils/sortd.c (revision d6ae5217716d0e83b63ef2baec5b10fcfb1fd4e8)
1e5c89e4eSSatish Balay /*
2e5c89e4eSSatish Balay    This file contains routines for sorting doubles.  Values are sorted in place.
3e5c89e4eSSatish Balay    These are provided because the general sort routines incur a great deal
4a5b23f4aSJose E. Roman    of overhead in calling the comparison routines.
5e5c89e4eSSatish Balay 
6e5c89e4eSSatish Balay  */
7c6db04a5SJed Brown #include <petscsys.h> /*I  "petscsys.h"  I*/
839f41f7fSStefano Zampini #include <petsc/private/petscimpl.h>
9e5c89e4eSSatish Balay 
109371c9d4SSatish Balay #define SWAP(a, b, t) \
11a8f51744SPierre Jolivet   do { \
129371c9d4SSatish Balay     t = a; \
139371c9d4SSatish Balay     a = b; \
149371c9d4SSatish Balay     b = t; \
15a8f51744SPierre Jolivet   } while (0)
16e5c89e4eSSatish Balay 
176a7c706bSVaclav Hapla /*@
18811af0c4SBarry Smith   PetscSortedReal - Determines whether the array of `PetscReal` is sorted.
196a7c706bSVaclav Hapla 
206a7c706bSVaclav Hapla   Not Collective
216a7c706bSVaclav Hapla 
226a7c706bSVaclav Hapla   Input Parameters:
236a7c706bSVaclav Hapla + n - number of values
246a7c706bSVaclav Hapla - X - array of integers
256a7c706bSVaclav Hapla 
262fe279fdSBarry Smith   Output Parameter:
276a7c706bSVaclav Hapla . sorted - flag whether the array is sorted
286a7c706bSVaclav Hapla 
296a7c706bSVaclav Hapla   Level: intermediate
306a7c706bSVaclav Hapla 
31db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortedInt()`, `PetscSortedMPIInt()`
326a7c706bSVaclav Hapla @*/
336497c311SBarry Smith PetscErrorCode PetscSortedReal(PetscCount n, const PetscReal X[], PetscBool *sorted)
34d71ae5a4SJacob Faibussowitsch {
356a7c706bSVaclav Hapla   PetscFunctionBegin;
366a7c706bSVaclav Hapla   PetscSorted(n, X, *sorted);
373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
386a7c706bSVaclav Hapla }
396a7c706bSVaclav Hapla 
40e5c89e4eSSatish Balay /* A simple version of quicksort; taken from Kernighan and Ritchie, page 87 */
416497c311SBarry Smith static PetscErrorCode PetscSortReal_Private(PetscReal *v, PetscCount right)
42d71ae5a4SJacob Faibussowitsch {
436497c311SBarry Smith   PetscCount i, last;
44e5c89e4eSSatish Balay   PetscReal  vl, tmp;
45e5c89e4eSSatish Balay 
46e5c89e4eSSatish Balay   PetscFunctionBegin;
47e5c89e4eSSatish Balay   if (right <= 1) {
48e5c89e4eSSatish Balay     if (right == 1) {
49e5c89e4eSSatish Balay       if (v[0] > v[1]) SWAP(v[0], v[1], tmp);
50e5c89e4eSSatish Balay     }
513ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
52e5c89e4eSSatish Balay   }
53e5c89e4eSSatish Balay   SWAP(v[0], v[right / 2], tmp);
54e5c89e4eSSatish Balay   vl   = v[0];
55e5c89e4eSSatish Balay   last = 0;
56e5c89e4eSSatish Balay   for (i = 1; i <= right; i++) {
579371c9d4SSatish Balay     if (v[i] < vl) {
589371c9d4SSatish Balay       last++;
599371c9d4SSatish Balay       SWAP(v[last], v[i], tmp);
609371c9d4SSatish Balay     }
61e5c89e4eSSatish Balay   }
62e5c89e4eSSatish Balay   SWAP(v[0], v[last], tmp);
633ba16761SJacob Faibussowitsch   PetscCall(PetscSortReal_Private(v, last - 1));
643ba16761SJacob Faibussowitsch   PetscCall(PetscSortReal_Private(v + last + 1, right - (last + 1)));
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66e5c89e4eSSatish Balay }
67e5c89e4eSSatish Balay 
68e5c89e4eSSatish Balay /*@
69811af0c4SBarry Smith   PetscSortReal - Sorts an array of `PetscReal` in place in increasing order.
70e5c89e4eSSatish Balay 
71e5c89e4eSSatish Balay   Not Collective
72e5c89e4eSSatish Balay 
73e5c89e4eSSatish Balay   Input Parameters:
74e5c89e4eSSatish Balay + n - number of values
75e5c89e4eSSatish Balay - v - array of doubles
76e5c89e4eSSatish Balay 
77667f096bSBarry Smith   Level: intermediate
78667f096bSBarry Smith 
79811af0c4SBarry Smith   Note:
80811af0c4SBarry Smith   This function serves as an alternative to `PetscRealSortSemiOrdered()`, and may perform faster especially if the array
81a5b23f4aSJose E. Roman   is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
82676f2a66SJacob Faibussowitsch   code to see which routine is fastest.
83676f2a66SJacob Faibussowitsch 
84db781477SPatrick Sanan .seealso: `PetscRealSortSemiOrdered()`, `PetscSortInt()`, `PetscSortRealWithPermutation()`, `PetscSortRealWithArrayInt()`
85e5c89e4eSSatish Balay @*/
866497c311SBarry Smith PetscErrorCode PetscSortReal(PetscCount n, PetscReal v[])
87d71ae5a4SJacob Faibussowitsch {
88e5c89e4eSSatish Balay   PetscFunctionBegin;
894f572ea9SToby Isaac   PetscAssertPointer(v, 2);
90e5c89e4eSSatish Balay   if (n < 8) {
916497c311SBarry Smith     PetscReal tmp, vk;
926497c311SBarry Smith     for (PetscCount k = 0; k < n; k++) {
93e5c89e4eSSatish Balay       vk = v[k];
946497c311SBarry Smith       for (PetscCount j = k + 1; j < n; j++) {
95e5c89e4eSSatish Balay         if (vk > v[j]) {
96e5c89e4eSSatish Balay           SWAP(v[k], v[j], tmp);
97e5c89e4eSSatish Balay           vk = v[k];
98e5c89e4eSSatish Balay         }
99e5c89e4eSSatish Balay       }
100e5c89e4eSSatish Balay     }
1016497c311SBarry Smith   } else {
1026497c311SBarry Smith     PetscInt N;
1036497c311SBarry Smith 
1046497c311SBarry Smith     PetscCall(PetscIntCast(n, &N));
1056497c311SBarry Smith     PetscCall(PetscSortReal_Private(v, N - 1));
1066497c311SBarry Smith   }
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108e5c89e4eSSatish Balay }
109e5c89e4eSSatish Balay 
1109371c9d4SSatish Balay #define SWAP2ri(a, b, c, d, rt, it) \
111a8f51744SPierre Jolivet   do { \
1129371c9d4SSatish Balay     rt = a; \
1139371c9d4SSatish Balay     a  = b; \
1149371c9d4SSatish Balay     b  = rt; \
1159371c9d4SSatish Balay     it = c; \
1169371c9d4SSatish Balay     c  = d; \
1179371c9d4SSatish Balay     d  = it; \
118a8f51744SPierre Jolivet   } while (0)
11939f41f7fSStefano Zampini 
12039f41f7fSStefano Zampini /* modified from PetscSortIntWithArray_Private */
1216497c311SBarry Smith static PetscErrorCode PetscSortRealWithArrayInt_Private(PetscReal *v, PetscInt *V, PetscCount right)
122d71ae5a4SJacob Faibussowitsch {
1236497c311SBarry Smith   PetscCount i, last;
1246497c311SBarry Smith   PetscInt   itmp;
12539f41f7fSStefano Zampini   PetscReal  rvl, rtmp;
12639f41f7fSStefano Zampini 
12739f41f7fSStefano Zampini   PetscFunctionBegin;
12839f41f7fSStefano Zampini   if (right <= 1) {
12939f41f7fSStefano Zampini     if (right == 1) {
13039f41f7fSStefano Zampini       if (v[0] > v[1]) SWAP2ri(v[0], v[1], V[0], V[1], rtmp, itmp);
13139f41f7fSStefano Zampini     }
1323ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
13339f41f7fSStefano Zampini   }
13439f41f7fSStefano Zampini   SWAP2ri(v[0], v[right / 2], V[0], V[right / 2], rtmp, itmp);
13539f41f7fSStefano Zampini   rvl  = v[0];
13639f41f7fSStefano Zampini   last = 0;
13739f41f7fSStefano Zampini   for (i = 1; i <= right; i++) {
1389371c9d4SSatish Balay     if (v[i] < rvl) {
1399371c9d4SSatish Balay       last++;
1409371c9d4SSatish Balay       SWAP2ri(v[last], v[i], V[last], V[i], rtmp, itmp);
1419371c9d4SSatish Balay     }
14239f41f7fSStefano Zampini   }
14339f41f7fSStefano Zampini   SWAP2ri(v[0], v[last], V[0], V[last], rtmp, itmp);
1449566063dSJacob Faibussowitsch   PetscCall(PetscSortRealWithArrayInt_Private(v, V, last - 1));
1459566063dSJacob Faibussowitsch   PetscCall(PetscSortRealWithArrayInt_Private(v + last + 1, V + last + 1, right - (last + 1)));
1463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14739f41f7fSStefano Zampini }
14839f41f7fSStefano Zampini /*@
149811af0c4SBarry Smith   PetscSortRealWithArrayInt - Sorts an array of `PetscReal` in place in increasing order;
150811af0c4SBarry Smith   changes a second `PetscInt` array to match the sorted first array.
15139f41f7fSStefano Zampini 
15239f41f7fSStefano Zampini   Not Collective
15339f41f7fSStefano Zampini 
15439f41f7fSStefano Zampini   Input Parameters:
15539f41f7fSStefano Zampini + n  - number of values
156aec76313SJacob Faibussowitsch . Ii - array of integers
157aec76313SJacob Faibussowitsch - r  - second array of integers
15839f41f7fSStefano Zampini 
15939f41f7fSStefano Zampini   Level: intermediate
16039f41f7fSStefano Zampini 
161db781477SPatrick Sanan .seealso: `PetscSortReal()`
16239f41f7fSStefano Zampini @*/
1636497c311SBarry Smith PetscErrorCode PetscSortRealWithArrayInt(PetscCount n, PetscReal r[], PetscInt Ii[])
164d71ae5a4SJacob Faibussowitsch {
1656497c311SBarry Smith   PetscCount j, k;
1666497c311SBarry Smith   PetscInt   itmp;
16739f41f7fSStefano Zampini   PetscReal  rk, rtmp;
16839f41f7fSStefano Zampini 
16939f41f7fSStefano Zampini   PetscFunctionBegin;
1704f572ea9SToby Isaac   PetscAssertPointer(r, 2);
1714f572ea9SToby Isaac   PetscAssertPointer(Ii, 3);
17239f41f7fSStefano Zampini   if (n < 8) {
17339f41f7fSStefano Zampini     for (k = 0; k < n; k++) {
17439f41f7fSStefano Zampini       rk = r[k];
17539f41f7fSStefano Zampini       for (j = k + 1; j < n; j++) {
17639f41f7fSStefano Zampini         if (rk > r[j]) {
17739f41f7fSStefano Zampini           SWAP2ri(r[k], r[j], Ii[k], Ii[j], rtmp, itmp);
17839f41f7fSStefano Zampini           rk = r[k];
17939f41f7fSStefano Zampini         }
18039f41f7fSStefano Zampini       }
18139f41f7fSStefano Zampini     }
18239f41f7fSStefano Zampini   } else {
1839566063dSJacob Faibussowitsch     PetscCall(PetscSortRealWithArrayInt_Private(r, Ii, n - 1));
18439f41f7fSStefano Zampini   }
1853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18639f41f7fSStefano Zampini }
18739f41f7fSStefano Zampini 
18839f41f7fSStefano Zampini /*@
1897a43aa34SPierre Jolivet   PetscFindReal - Finds a `PetscReal` in a sorted array of `PetscReal`s
19039f41f7fSStefano Zampini 
19139f41f7fSStefano Zampini   Not Collective
19239f41f7fSStefano Zampini 
19339f41f7fSStefano Zampini   Input Parameters:
19439f41f7fSStefano Zampini + key - the value to locate
19539f41f7fSStefano Zampini . n   - number of values in the array
196aec76313SJacob Faibussowitsch . t   - array of values
19739f41f7fSStefano Zampini - eps - tolerance used to compare
19839f41f7fSStefano Zampini 
19939f41f7fSStefano Zampini   Output Parameter:
20039f41f7fSStefano Zampini . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
20139f41f7fSStefano Zampini 
20239f41f7fSStefano Zampini   Level: intermediate
20339f41f7fSStefano Zampini 
204db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortRealWithArrayInt()`
20539f41f7fSStefano Zampini @*/
2066497c311SBarry Smith PetscErrorCode PetscFindReal(PetscReal key, PetscCount n, const PetscReal t[], PetscReal eps, PetscInt *loc)
207d71ae5a4SJacob Faibussowitsch {
2086497c311SBarry Smith   PetscInt lo = 0, hi;
20939f41f7fSStefano Zampini 
21039f41f7fSStefano Zampini   PetscFunctionBegin;
2114f572ea9SToby Isaac   PetscAssertPointer(loc, 5);
2129371c9d4SSatish Balay   if (!n) {
2139371c9d4SSatish Balay     *loc = -1;
2143ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2159371c9d4SSatish Balay   }
2164f572ea9SToby Isaac   PetscAssertPointer(t, 3);
2176497c311SBarry Smith   PetscCall(PetscIntCast(n, &hi));
21839f41f7fSStefano Zampini   while (hi - lo > 1) {
21939f41f7fSStefano Zampini     PetscInt mid = lo + (hi - lo) / 2;
22042f965a0SMatthew G. Knepley     PetscAssert(t[lo] <= t[mid] && t[mid] <= t[hi - 1], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input array was not sorted: (%g, %g, %g)", (double)t[lo], (double)t[mid], (double)t[hi - 1]);
22139f41f7fSStefano Zampini     if (key < t[mid]) hi = mid;
22239f41f7fSStefano Zampini     else lo = mid;
22339f41f7fSStefano Zampini   }
22439f41f7fSStefano Zampini   *loc = (PetscAbsReal(key - t[lo]) < eps) ? lo : -(lo + (key > t[lo]) + 1);
2253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22639f41f7fSStefano Zampini }
227745b41b2SMatthew G. Knepley 
228745b41b2SMatthew G. Knepley /*@
229811af0c4SBarry Smith   PetscSortRemoveDupsReal - Sorts an array of `PetscReal` in place in increasing order and removes all duplicate entries
230745b41b2SMatthew G. Knepley 
231745b41b2SMatthew G. Knepley   Not Collective
232745b41b2SMatthew G. Knepley 
233745b41b2SMatthew G. Knepley   Input Parameters:
23426a11704SBarry Smith + n - initial number of values
23526a11704SBarry Smith - v - array of values
236745b41b2SMatthew G. Knepley 
237*d6ae5217SJose E. Roman   Note:
238*d6ae5217SJose E. Roman   On output both `n` and `v` are modified with non-redundant values.
239745b41b2SMatthew G. Knepley 
240745b41b2SMatthew G. Knepley   Level: intermediate
241745b41b2SMatthew G. Knepley 
242db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortRemoveDupsInt()`
243745b41b2SMatthew G. Knepley @*/
244d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRemoveDupsReal(PetscInt *n, PetscReal v[])
245d71ae5a4SJacob Faibussowitsch {
246745b41b2SMatthew G. Knepley   PetscInt i, s = 0, N = *n, b = 0;
247745b41b2SMatthew G. Knepley 
248745b41b2SMatthew G. Knepley   PetscFunctionBegin;
2499566063dSJacob Faibussowitsch   PetscCall(PetscSortReal(N, v));
250745b41b2SMatthew G. Knepley   for (i = 0; i < N - 1; i++) {
251745b41b2SMatthew G. Knepley     if (v[b + s + 1] != v[b]) {
2529371c9d4SSatish Balay       v[b + 1] = v[b + s + 1];
2539371c9d4SSatish Balay       b++;
254745b41b2SMatthew G. Knepley     } else s++;
255745b41b2SMatthew G. Knepley   }
256745b41b2SMatthew G. Knepley   *n = N - s;
2573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
258745b41b2SMatthew G. Knepley }
259745b41b2SMatthew G. Knepley 
260d97c5584SHong Zhang /*@
261811af0c4SBarry Smith   PetscSortSplit - Quick-sort split of an array of `PetscScalar`s in place.
262d97c5584SHong Zhang 
263d97c5584SHong Zhang   Not Collective
264d97c5584SHong Zhang 
265d97c5584SHong Zhang   Input Parameters:
2667a43aa34SPierre Jolivet + ncut - splitting index
2676b867d5aSJose E. Roman - n    - number of values to sort
268d97c5584SHong Zhang 
2696b867d5aSJose E. Roman   Input/Output Parameters:
2706b867d5aSJose E. Roman + a   - array of values, on output the values are permuted such that its elements satisfy:
271d97c5584SHong Zhang            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
272d97c5584SHong Zhang            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
2736b867d5aSJose E. Roman - idx - index for array a, on output permuted accordingly
274d97c5584SHong Zhang 
275d97c5584SHong Zhang   Level: intermediate
276d97c5584SHong Zhang 
277db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()`
278d97c5584SHong Zhang @*/
279d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortSplit(PetscInt ncut, PetscInt n, PetscScalar a[], PetscInt idx[])
280d71ae5a4SJacob Faibussowitsch {
281d97c5584SHong Zhang   PetscInt    i, mid, last, itmp, j, first;
282d97c5584SHong Zhang   PetscScalar d, tmp;
283d97c5584SHong Zhang   PetscReal   abskey;
284d97c5584SHong Zhang 
285d97c5584SHong Zhang   PetscFunctionBegin;
286d97c5584SHong Zhang   first = 0;
287d97c5584SHong Zhang   last  = n - 1;
2883ba16761SJacob Faibussowitsch   if (ncut < first || ncut > last) PetscFunctionReturn(PETSC_SUCCESS);
289d97c5584SHong Zhang 
290d97c5584SHong Zhang   while (1) {
291d97c5584SHong Zhang     mid    = first;
2922a274a78SSatish Balay     d      = a[mid];
2932a274a78SSatish Balay     abskey = PetscAbsScalar(d);
294d97c5584SHong Zhang     i      = last;
295d97c5584SHong Zhang     for (j = first + 1; j <= i; ++j) {
2962a274a78SSatish Balay       d = a[j];
2972a274a78SSatish Balay       if (PetscAbsScalar(d) >= abskey) {
298d97c5584SHong Zhang         ++mid;
299d97c5584SHong Zhang         /* interchange */
3009371c9d4SSatish Balay         tmp      = a[mid];
3019371c9d4SSatish Balay         itmp     = idx[mid];
3029371c9d4SSatish Balay         a[mid]   = a[j];
3039371c9d4SSatish Balay         idx[mid] = idx[j];
3049371c9d4SSatish Balay         a[j]     = tmp;
3059371c9d4SSatish Balay         idx[j]   = itmp;
306d97c5584SHong Zhang       }
307d97c5584SHong Zhang     }
308d97c5584SHong Zhang 
309d97c5584SHong Zhang     /* interchange */
3109371c9d4SSatish Balay     tmp        = a[mid];
3119371c9d4SSatish Balay     itmp       = idx[mid];
3129371c9d4SSatish Balay     a[mid]     = a[first];
3139371c9d4SSatish Balay     idx[mid]   = idx[first];
3149371c9d4SSatish Balay     a[first]   = tmp;
3159371c9d4SSatish Balay     idx[first] = itmp;
316d97c5584SHong Zhang 
317d97c5584SHong Zhang     /* test for while loop */
318a297a907SKarl Rupp     if (mid == ncut) break;
319a297a907SKarl Rupp     else if (mid > ncut) last = mid - 1;
320a297a907SKarl Rupp     else first = mid + 1;
321d97c5584SHong Zhang   }
3223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
323d97c5584SHong Zhang }
324021822bcSHong Zhang 
325021822bcSHong Zhang /*@
326811af0c4SBarry Smith   PetscSortSplitReal - Quick-sort split of an array of `PetscReal`s in place.
327021822bcSHong Zhang 
328021822bcSHong Zhang   Not Collective
329021822bcSHong Zhang 
330021822bcSHong Zhang   Input Parameters:
3317a43aa34SPierre Jolivet + ncut - splitting index
3326b867d5aSJose E. Roman - n    - number of values to sort
333021822bcSHong Zhang 
3346b867d5aSJose E. Roman   Input/Output Parameters:
3356b867d5aSJose E. Roman + a   - array of values, on output the values are permuted such that its elements satisfy:
336021822bcSHong Zhang            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
337021822bcSHong Zhang            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
3386b867d5aSJose E. Roman - idx - index for array a, on output permuted accordingly
339021822bcSHong Zhang 
340021822bcSHong Zhang   Level: intermediate
341021822bcSHong Zhang 
342db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()`
343021822bcSHong Zhang @*/
344d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortSplitReal(PetscInt ncut, PetscInt n, PetscReal a[], PetscInt idx[])
345d71ae5a4SJacob Faibussowitsch {
346021822bcSHong Zhang   PetscInt  i, mid, last, itmp, j, first;
347021822bcSHong Zhang   PetscReal d, tmp;
348021822bcSHong Zhang   PetscReal abskey;
349021822bcSHong Zhang 
350021822bcSHong Zhang   PetscFunctionBegin;
351021822bcSHong Zhang   first = 0;
352021822bcSHong Zhang   last  = n - 1;
3533ba16761SJacob Faibussowitsch   if (ncut < first || ncut > last) PetscFunctionReturn(PETSC_SUCCESS);
354021822bcSHong Zhang 
355021822bcSHong Zhang   while (1) {
356021822bcSHong Zhang     mid    = first;
3572a274a78SSatish Balay     d      = a[mid];
3582a274a78SSatish Balay     abskey = PetscAbsReal(d);
359021822bcSHong Zhang     i      = last;
360021822bcSHong Zhang     for (j = first + 1; j <= i; ++j) {
3612a274a78SSatish Balay       d = a[j];
3622a274a78SSatish Balay       if (PetscAbsReal(d) >= abskey) {
363021822bcSHong Zhang         ++mid;
364021822bcSHong Zhang         /* interchange */
3659371c9d4SSatish Balay         tmp      = a[mid];
3669371c9d4SSatish Balay         itmp     = idx[mid];
3679371c9d4SSatish Balay         a[mid]   = a[j];
3689371c9d4SSatish Balay         idx[mid] = idx[j];
3699371c9d4SSatish Balay         a[j]     = tmp;
3709371c9d4SSatish Balay         idx[j]   = itmp;
371021822bcSHong Zhang       }
372021822bcSHong Zhang     }
373021822bcSHong Zhang 
374021822bcSHong Zhang     /* interchange */
3759371c9d4SSatish Balay     tmp        = a[mid];
3769371c9d4SSatish Balay     itmp       = idx[mid];
3779371c9d4SSatish Balay     a[mid]     = a[first];
3789371c9d4SSatish Balay     idx[mid]   = idx[first];
3799371c9d4SSatish Balay     a[first]   = tmp;
3809371c9d4SSatish Balay     idx[first] = itmp;
381021822bcSHong Zhang 
382021822bcSHong Zhang     /* test for while loop */
383a297a907SKarl Rupp     if (mid == ncut) break;
384a297a907SKarl Rupp     else if (mid > ncut) last = mid - 1;
385a297a907SKarl Rupp     else first = mid + 1;
386021822bcSHong Zhang   }
3873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
388021822bcSHong Zhang }
389