xref: /petsc/src/sys/utils/sortd.c (revision 676f2a665bc48dc8c31863fa9f09556a43eff7ab)
17d0a6c19SBarry Smith 
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file contains routines for sorting doubles.  Values are sorted in place.
4e5c89e4eSSatish Balay    These are provided because the general sort routines incur a great deal
5e5c89e4eSSatish Balay    of overhead in calling the comparision routines.
6e5c89e4eSSatish Balay 
7e5c89e4eSSatish Balay  */
8c6db04a5SJed Brown #include <petscsys.h>                /*I  "petscsys.h"  I*/
939f41f7fSStefano Zampini #include <petsc/private/petscimpl.h>
10e5c89e4eSSatish Balay 
11e5c89e4eSSatish Balay #define SWAP(a,b,t) {t=a;a=b;b=t;}
12e5c89e4eSSatish Balay 
136a7c706bSVaclav Hapla /*@
146a7c706bSVaclav Hapla    PetscSortedReal - Determines whether the array is sorted.
156a7c706bSVaclav Hapla 
166a7c706bSVaclav Hapla    Not Collective
176a7c706bSVaclav Hapla 
186a7c706bSVaclav Hapla    Input Parameters:
196a7c706bSVaclav Hapla +  n  - number of values
206a7c706bSVaclav Hapla -  X  - array of integers
216a7c706bSVaclav Hapla 
226a7c706bSVaclav Hapla    Output Parameters:
236a7c706bSVaclav Hapla .  sorted - flag whether the array is sorted
246a7c706bSVaclav Hapla 
256a7c706bSVaclav Hapla    Level: intermediate
266a7c706bSVaclav Hapla 
276a7c706bSVaclav Hapla .seealso: PetscSortReal(), PetscSortedInt(), PetscSortedMPIInt()
286a7c706bSVaclav Hapla @*/
296a7c706bSVaclav Hapla PetscErrorCode  PetscSortedReal(PetscInt n,const PetscReal X[],PetscBool *sorted)
306a7c706bSVaclav Hapla {
316a7c706bSVaclav Hapla   PetscFunctionBegin;
326a7c706bSVaclav Hapla   PetscSorted(n,X,*sorted);
336a7c706bSVaclav Hapla   PetscFunctionReturn(0);
346a7c706bSVaclav Hapla }
356a7c706bSVaclav Hapla 
36e5c89e4eSSatish Balay /* A simple version of quicksort; taken from Kernighan and Ritchie, page 87 */
37e5c89e4eSSatish Balay static PetscErrorCode PetscSortReal_Private(PetscReal *v,PetscInt right)
38e5c89e4eSSatish Balay {
39e5c89e4eSSatish Balay   PetscInt  i,last;
40e5c89e4eSSatish Balay   PetscReal vl,tmp;
41e5c89e4eSSatish Balay 
42e5c89e4eSSatish Balay   PetscFunctionBegin;
43e5c89e4eSSatish Balay   if (right <= 1) {
44e5c89e4eSSatish Balay     if (right == 1) {
45e5c89e4eSSatish Balay       if (v[0] > v[1]) SWAP(v[0],v[1],tmp);
46e5c89e4eSSatish Balay     }
47e5c89e4eSSatish Balay     PetscFunctionReturn(0);
48e5c89e4eSSatish Balay   }
49e5c89e4eSSatish Balay   SWAP(v[0],v[right/2],tmp);
50e5c89e4eSSatish Balay   vl   = v[0];
51e5c89e4eSSatish Balay   last = 0;
52e5c89e4eSSatish Balay   for (i=1; i<=right; i++) {
53e5c89e4eSSatish Balay     if (v[i] < vl) {last++; SWAP(v[last],v[i],tmp);}
54e5c89e4eSSatish Balay   }
55e5c89e4eSSatish Balay   SWAP(v[0],v[last],tmp);
56e5c89e4eSSatish Balay   PetscSortReal_Private(v,last-1);
57e5c89e4eSSatish Balay   PetscSortReal_Private(v+last+1,right-(last+1));
58e5c89e4eSSatish Balay   PetscFunctionReturn(0);
59e5c89e4eSSatish Balay }
60e5c89e4eSSatish Balay 
61e5c89e4eSSatish Balay /*@
62e5c89e4eSSatish Balay    PetscSortReal - Sorts an array of doubles in place in increasing order.
63e5c89e4eSSatish Balay 
64e5c89e4eSSatish Balay    Not Collective
65e5c89e4eSSatish Balay 
66e5c89e4eSSatish Balay    Input Parameters:
67e5c89e4eSSatish Balay +  n  - number of values
68e5c89e4eSSatish Balay -  v  - array of doubles
69e5c89e4eSSatish Balay 
70*676f2a66SJacob Faibussowitsch    Notes:
71*676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscRealSortSemiOrdered(), and may perform faster especially if the array
72*676f2a66SJacob Faibussowitsch    is completely random. There are exceptions to this and so it is __highly__ recomended that the user benchmark their
73*676f2a66SJacob Faibussowitsch    code to see which routine is fastest.
74*676f2a66SJacob Faibussowitsch 
75e5c89e4eSSatish Balay    Level: intermediate
76e5c89e4eSSatish Balay 
77*676f2a66SJacob Faibussowitsch .seealso: PetscRealSortSemiOrdered(), PetscSortInt(), PetscSortRealWithPermutation(), PetscSortRealWithArrayInt()
78e5c89e4eSSatish Balay @*/
797087cfbeSBarry Smith PetscErrorCode  PetscSortReal(PetscInt n,PetscReal v[])
80e5c89e4eSSatish Balay {
81e5c89e4eSSatish Balay   PetscInt  j,k;
82e5c89e4eSSatish Balay   PetscReal tmp,vk;
83e5c89e4eSSatish Balay 
84e5c89e4eSSatish Balay   PetscFunctionBegin;
8539f41f7fSStefano Zampini   PetscValidPointer(v,2);
86e5c89e4eSSatish Balay   if (n<8) {
87e5c89e4eSSatish Balay     for (k=0; k<n; k++) {
88e5c89e4eSSatish Balay       vk = v[k];
89e5c89e4eSSatish Balay       for (j=k+1; j<n; j++) {
90e5c89e4eSSatish Balay         if (vk > v[j]) {
91e5c89e4eSSatish Balay           SWAP(v[k],v[j],tmp);
92e5c89e4eSSatish Balay           vk = v[k];
93e5c89e4eSSatish Balay         }
94e5c89e4eSSatish Balay       }
95e5c89e4eSSatish Balay     }
96a297a907SKarl Rupp   } else PetscSortReal_Private(v,n-1);
97e5c89e4eSSatish Balay   PetscFunctionReturn(0);
98e5c89e4eSSatish Balay }
99e5c89e4eSSatish Balay 
10039f41f7fSStefano Zampini #define SWAP2ri(a,b,c,d,rt,it) {rt=a;a=b;b=rt;it=c;c=d;d=it;}
10139f41f7fSStefano Zampini 
10239f41f7fSStefano Zampini /* modified from PetscSortIntWithArray_Private */
10339f41f7fSStefano Zampini static PetscErrorCode PetscSortRealWithArrayInt_Private(PetscReal *v,PetscInt *V,PetscInt right)
10439f41f7fSStefano Zampini {
10539f41f7fSStefano Zampini   PetscErrorCode ierr;
10639f41f7fSStefano Zampini   PetscInt       i,last,itmp;
10739f41f7fSStefano Zampini   PetscReal      rvl,rtmp;
10839f41f7fSStefano Zampini 
10939f41f7fSStefano Zampini   PetscFunctionBegin;
11039f41f7fSStefano Zampini   if (right <= 1) {
11139f41f7fSStefano Zampini     if (right == 1) {
11239f41f7fSStefano Zampini       if (v[0] > v[1]) SWAP2ri(v[0],v[1],V[0],V[1],rtmp,itmp);
11339f41f7fSStefano Zampini     }
11439f41f7fSStefano Zampini     PetscFunctionReturn(0);
11539f41f7fSStefano Zampini   }
11639f41f7fSStefano Zampini   SWAP2ri(v[0],v[right/2],V[0],V[right/2],rtmp,itmp);
11739f41f7fSStefano Zampini   rvl  = v[0];
11839f41f7fSStefano Zampini   last = 0;
11939f41f7fSStefano Zampini   for (i=1; i<=right; i++) {
12039f41f7fSStefano Zampini     if (v[i] < rvl) {last++; SWAP2ri(v[last],v[i],V[last],V[i],rtmp,itmp);}
12139f41f7fSStefano Zampini   }
12239f41f7fSStefano Zampini   SWAP2ri(v[0],v[last],V[0],V[last],rtmp,itmp);
12339f41f7fSStefano Zampini   ierr = PetscSortRealWithArrayInt_Private(v,V,last-1);CHKERRQ(ierr);
12439f41f7fSStefano Zampini   ierr = PetscSortRealWithArrayInt_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr);
12539f41f7fSStefano Zampini   PetscFunctionReturn(0);
12639f41f7fSStefano Zampini }
12739f41f7fSStefano Zampini /*@
12839f41f7fSStefano Zampini    PetscSortRealWithArrayInt - Sorts an array of PetscReal in place in increasing order;
12939f41f7fSStefano Zampini        changes a second PetscInt array to match the sorted first array.
13039f41f7fSStefano Zampini 
13139f41f7fSStefano Zampini    Not Collective
13239f41f7fSStefano Zampini 
13339f41f7fSStefano Zampini    Input Parameters:
13439f41f7fSStefano Zampini +  n  - number of values
13539f41f7fSStefano Zampini .  i  - array of integers
13639f41f7fSStefano Zampini -  I - second array of integers
13739f41f7fSStefano Zampini 
13839f41f7fSStefano Zampini    Level: intermediate
13939f41f7fSStefano Zampini 
14039f41f7fSStefano Zampini .seealso: PetscSortReal()
14139f41f7fSStefano Zampini @*/
14239f41f7fSStefano Zampini PetscErrorCode  PetscSortRealWithArrayInt(PetscInt n,PetscReal r[],PetscInt Ii[])
14339f41f7fSStefano Zampini {
14439f41f7fSStefano Zampini   PetscErrorCode ierr;
14539f41f7fSStefano Zampini   PetscInt       j,k,itmp;
14639f41f7fSStefano Zampini   PetscReal      rk,rtmp;
14739f41f7fSStefano Zampini 
14839f41f7fSStefano Zampini   PetscFunctionBegin;
14939f41f7fSStefano Zampini   PetscValidPointer(r,2);
15039f41f7fSStefano Zampini   PetscValidPointer(Ii,3);
15139f41f7fSStefano Zampini   if (n<8) {
15239f41f7fSStefano Zampini     for (k=0; k<n; k++) {
15339f41f7fSStefano Zampini       rk = r[k];
15439f41f7fSStefano Zampini       for (j=k+1; j<n; j++) {
15539f41f7fSStefano Zampini         if (rk > r[j]) {
15639f41f7fSStefano Zampini           SWAP2ri(r[k],r[j],Ii[k],Ii[j],rtmp,itmp);
15739f41f7fSStefano Zampini           rk = r[k];
15839f41f7fSStefano Zampini         }
15939f41f7fSStefano Zampini       }
16039f41f7fSStefano Zampini     }
16139f41f7fSStefano Zampini   } else {
16239f41f7fSStefano Zampini     ierr = PetscSortRealWithArrayInt_Private(r,Ii,n-1);CHKERRQ(ierr);
16339f41f7fSStefano Zampini   }
16439f41f7fSStefano Zampini   PetscFunctionReturn(0);
16539f41f7fSStefano Zampini }
16639f41f7fSStefano Zampini 
16739f41f7fSStefano Zampini /*@
16839f41f7fSStefano Zampini   PetscFindReal - Finds a PetscReal in a sorted array of PetscReals
16939f41f7fSStefano Zampini 
17039f41f7fSStefano Zampini    Not Collective
17139f41f7fSStefano Zampini 
17239f41f7fSStefano Zampini    Input Parameters:
17339f41f7fSStefano Zampini +  key - the value to locate
17439f41f7fSStefano Zampini .  n   - number of values in the array
17539f41f7fSStefano Zampini .  ii  - array of values
17639f41f7fSStefano Zampini -  eps - tolerance used to compare
17739f41f7fSStefano Zampini 
17839f41f7fSStefano Zampini    Output Parameter:
17939f41f7fSStefano Zampini .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
18039f41f7fSStefano Zampini 
18139f41f7fSStefano Zampini    Level: intermediate
18239f41f7fSStefano Zampini 
18339f41f7fSStefano Zampini .seealso: PetscSortReal(), PetscSortRealWithArrayInt()
18439f41f7fSStefano Zampini @*/
18539f41f7fSStefano Zampini PetscErrorCode PetscFindReal(PetscReal key, PetscInt n, const PetscReal t[], PetscReal eps, PetscInt *loc)
18639f41f7fSStefano Zampini {
18739f41f7fSStefano Zampini   PetscInt lo = 0,hi = n;
18839f41f7fSStefano Zampini 
18939f41f7fSStefano Zampini   PetscFunctionBegin;
19039f41f7fSStefano Zampini   PetscValidPointer(loc,4);
19139f41f7fSStefano Zampini   if (!n) {*loc = -1; PetscFunctionReturn(0);}
19239f41f7fSStefano Zampini   PetscValidPointer(t,3);
1936a7c706bSVaclav Hapla   PetscCheckSorted(n,t);
19439f41f7fSStefano Zampini   while (hi - lo > 1) {
19539f41f7fSStefano Zampini     PetscInt mid = lo + (hi - lo)/2;
19639f41f7fSStefano Zampini     if (key < t[mid]) hi = mid;
19739f41f7fSStefano Zampini     else              lo = mid;
19839f41f7fSStefano Zampini   }
19939f41f7fSStefano Zampini   *loc = (PetscAbsReal(key - t[lo]) < eps) ? lo : -(lo + (key > t[lo]) + 1);
20039f41f7fSStefano Zampini   PetscFunctionReturn(0);
20139f41f7fSStefano Zampini }
202745b41b2SMatthew G. Knepley 
203745b41b2SMatthew G. Knepley /*@
204745b41b2SMatthew G. Knepley    PetscSortRemoveDupsReal - Sorts an array of doubles in place in increasing order removes all duplicate entries
205745b41b2SMatthew G. Knepley 
206745b41b2SMatthew G. Knepley    Not Collective
207745b41b2SMatthew G. Knepley 
208745b41b2SMatthew G. Knepley    Input Parameters:
209745b41b2SMatthew G. Knepley +  n  - number of values
210745b41b2SMatthew G. Knepley -  v  - array of doubles
211745b41b2SMatthew G. Knepley 
212745b41b2SMatthew G. Knepley    Output Parameter:
213745b41b2SMatthew G. Knepley .  n - number of non-redundant values
214745b41b2SMatthew G. Knepley 
215745b41b2SMatthew G. Knepley    Level: intermediate
216745b41b2SMatthew G. Knepley 
217745b41b2SMatthew G. Knepley .seealso: PetscSortReal(), PetscSortRemoveDupsInt()
218745b41b2SMatthew G. Knepley @*/
219745b41b2SMatthew G. Knepley PetscErrorCode  PetscSortRemoveDupsReal(PetscInt *n,PetscReal v[])
220745b41b2SMatthew G. Knepley {
221745b41b2SMatthew G. Knepley   PetscErrorCode ierr;
222745b41b2SMatthew G. Knepley   PetscInt       i,s = 0,N = *n, b = 0;
223745b41b2SMatthew G. Knepley 
224745b41b2SMatthew G. Knepley   PetscFunctionBegin;
225745b41b2SMatthew G. Knepley   ierr = PetscSortReal(N,v);CHKERRQ(ierr);
226745b41b2SMatthew G. Knepley   for (i=0; i<N-1; i++) {
227745b41b2SMatthew G. Knepley     if (v[b+s+1] != v[b]) {
228745b41b2SMatthew G. Knepley       v[b+1] = v[b+s+1]; b++;
229745b41b2SMatthew G. Knepley     } else s++;
230745b41b2SMatthew G. Knepley   }
231745b41b2SMatthew G. Knepley   *n = N - s;
232745b41b2SMatthew G. Knepley   PetscFunctionReturn(0);
233745b41b2SMatthew G. Knepley }
234745b41b2SMatthew G. Knepley 
235d97c5584SHong Zhang /*@
236021822bcSHong Zhang    PetscSortSplit - Quick-sort split of an array of PetscScalars in place.
237d97c5584SHong Zhang 
238d97c5584SHong Zhang    Not Collective
239d97c5584SHong Zhang 
240d97c5584SHong Zhang    Input Parameters:
241d97c5584SHong Zhang +  ncut  - splitig index
242d97c5584SHong Zhang .  n     - number of values to sort
243d97c5584SHong Zhang .  a     - array of values
244d97c5584SHong Zhang -  idx   - index for array a
245d97c5584SHong Zhang 
246d97c5584SHong Zhang    Output Parameters:
247d97c5584SHong Zhang +  a     - permuted array of values such that its elements satisfy:
248d97c5584SHong Zhang            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
249d97c5584SHong Zhang            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
250d97c5584SHong Zhang -  idx   - permuted index of array a
251d97c5584SHong Zhang 
252d97c5584SHong Zhang    Level: intermediate
253d97c5584SHong Zhang 
254d97c5584SHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation()
255d97c5584SHong Zhang @*/
2567087cfbeSBarry Smith PetscErrorCode  PetscSortSplit(PetscInt ncut,PetscInt n,PetscScalar a[],PetscInt idx[])
257d97c5584SHong Zhang {
258d97c5584SHong Zhang   PetscInt    i,mid,last,itmp,j,first;
259d97c5584SHong Zhang   PetscScalar d,tmp;
260d97c5584SHong Zhang   PetscReal   abskey;
261d97c5584SHong Zhang 
262d97c5584SHong Zhang   PetscFunctionBegin;
263d97c5584SHong Zhang   first = 0;
264d97c5584SHong Zhang   last  = n-1;
265d97c5584SHong Zhang   if (ncut < first || ncut > last) PetscFunctionReturn(0);
266d97c5584SHong Zhang 
267d97c5584SHong Zhang   while (1) {
268d97c5584SHong Zhang     mid    = first;
2692a274a78SSatish Balay     d      = a[mid];
2702a274a78SSatish Balay     abskey = PetscAbsScalar(d);
271d97c5584SHong Zhang     i      = last;
272d97c5584SHong Zhang     for (j = first + 1; j <= i; ++j) {
2732a274a78SSatish Balay       d = a[j];
2742a274a78SSatish Balay       if (PetscAbsScalar(d) >= abskey) {
275d97c5584SHong Zhang         ++mid;
276d97c5584SHong Zhang         /* interchange */
277d97c5584SHong Zhang         tmp = a[mid];  itmp = idx[mid];
278d97c5584SHong Zhang         a[mid] = a[j]; idx[mid] = idx[j];
279d97c5584SHong Zhang         a[j] = tmp;    idx[j] = itmp;
280d97c5584SHong Zhang       }
281d97c5584SHong Zhang     }
282d97c5584SHong Zhang 
283d97c5584SHong Zhang     /* interchange */
284d97c5584SHong Zhang     tmp = a[mid];      itmp = idx[mid];
285d97c5584SHong Zhang     a[mid] = a[first]; idx[mid] = idx[first];
286d97c5584SHong Zhang     a[first] = tmp;    idx[first] = itmp;
287d97c5584SHong Zhang 
288d97c5584SHong Zhang     /* test for while loop */
289a297a907SKarl Rupp     if (mid == ncut) break;
290a297a907SKarl Rupp     else if (mid > ncut) last = mid - 1;
291a297a907SKarl Rupp     else first = mid + 1;
292d97c5584SHong Zhang   }
293d97c5584SHong Zhang   PetscFunctionReturn(0);
294d97c5584SHong Zhang }
295021822bcSHong Zhang 
296021822bcSHong Zhang /*@
297021822bcSHong Zhang    PetscSortSplitReal - Quick-sort split of an array of PetscReals in place.
298021822bcSHong Zhang 
299021822bcSHong Zhang    Not Collective
300021822bcSHong Zhang 
301021822bcSHong Zhang    Input Parameters:
302021822bcSHong Zhang +  ncut  - splitig index
303021822bcSHong Zhang .  n     - number of values to sort
304021822bcSHong Zhang .  a     - array of values in PetscReal
305021822bcSHong Zhang -  idx   - index for array a
306021822bcSHong Zhang 
307021822bcSHong Zhang    Output Parameters:
308021822bcSHong Zhang +  a     - permuted array of real values such that its elements satisfy:
309021822bcSHong Zhang            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
310021822bcSHong Zhang            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
311021822bcSHong Zhang -  idx   - permuted index of array a
312021822bcSHong Zhang 
313021822bcSHong Zhang    Level: intermediate
314021822bcSHong Zhang 
315021822bcSHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation()
316021822bcSHong Zhang @*/
3177087cfbeSBarry Smith PetscErrorCode  PetscSortSplitReal(PetscInt ncut,PetscInt n,PetscReal a[],PetscInt idx[])
318021822bcSHong Zhang {
319021822bcSHong Zhang   PetscInt  i,mid,last,itmp,j,first;
320021822bcSHong Zhang   PetscReal d,tmp;
321021822bcSHong Zhang   PetscReal abskey;
322021822bcSHong Zhang 
323021822bcSHong Zhang   PetscFunctionBegin;
324021822bcSHong Zhang   first = 0;
325021822bcSHong Zhang   last  = n-1;
326021822bcSHong Zhang   if (ncut < first || ncut > last) PetscFunctionReturn(0);
327021822bcSHong Zhang 
328021822bcSHong Zhang   while (1) {
329021822bcSHong Zhang     mid    = first;
3302a274a78SSatish Balay     d      = a[mid];
3312a274a78SSatish Balay     abskey = PetscAbsReal(d);
332021822bcSHong Zhang     i      = last;
333021822bcSHong Zhang     for (j = first + 1; j <= i; ++j) {
3342a274a78SSatish Balay       d = a[j];
3352a274a78SSatish Balay       if (PetscAbsReal(d) >= abskey) {
336021822bcSHong Zhang         ++mid;
337021822bcSHong Zhang         /* interchange */
338021822bcSHong Zhang         tmp = a[mid];  itmp = idx[mid];
339021822bcSHong Zhang         a[mid] = a[j]; idx[mid] = idx[j];
340021822bcSHong Zhang         a[j] = tmp;    idx[j] = itmp;
341021822bcSHong Zhang       }
342021822bcSHong Zhang     }
343021822bcSHong Zhang 
344021822bcSHong Zhang     /* interchange */
345021822bcSHong Zhang     tmp = a[mid];      itmp = idx[mid];
346021822bcSHong Zhang     a[mid] = a[first]; idx[mid] = idx[first];
347021822bcSHong Zhang     a[first] = tmp;    idx[first] = itmp;
348021822bcSHong Zhang 
349021822bcSHong Zhang     /* test for while loop */
350a297a907SKarl Rupp     if (mid == ncut) break;
351a297a907SKarl Rupp     else if (mid > ncut) last = mid - 1;
352a297a907SKarl Rupp     else first = mid + 1;
353021822bcSHong Zhang   }
354021822bcSHong Zhang   PetscFunctionReturn(0);
355021822bcSHong Zhang }
356