xref: /petsc/src/sys/utils/sortd.c (revision dadcf80911fb48939c55327431ae8d7e47dbe367)
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
5a5b23f4aSJose E. Roman    of overhead in calling the comparison 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 
70676f2a66SJacob Faibussowitsch    Notes:
71676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscRealSortSemiOrdered(), and may perform faster especially if the array
72a5b23f4aSJose E. Roman    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
73676f2a66SJacob Faibussowitsch    code to see which routine is fastest.
74676f2a66SJacob Faibussowitsch 
75e5c89e4eSSatish Balay    Level: intermediate
76e5c89e4eSSatish Balay 
77676f2a66SJacob 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;
85*dadcf809SJacob Faibussowitsch   PetscValidRealPointer(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   PetscInt       i,last,itmp;
10639f41f7fSStefano Zampini   PetscReal      rvl,rtmp;
10739f41f7fSStefano Zampini 
10839f41f7fSStefano Zampini   PetscFunctionBegin;
10939f41f7fSStefano Zampini   if (right <= 1) {
11039f41f7fSStefano Zampini     if (right == 1) {
11139f41f7fSStefano Zampini       if (v[0] > v[1]) SWAP2ri(v[0],v[1],V[0],V[1],rtmp,itmp);
11239f41f7fSStefano Zampini     }
11339f41f7fSStefano Zampini     PetscFunctionReturn(0);
11439f41f7fSStefano Zampini   }
11539f41f7fSStefano Zampini   SWAP2ri(v[0],v[right/2],V[0],V[right/2],rtmp,itmp);
11639f41f7fSStefano Zampini   rvl  = v[0];
11739f41f7fSStefano Zampini   last = 0;
11839f41f7fSStefano Zampini   for (i=1; i<=right; i++) {
11939f41f7fSStefano Zampini     if (v[i] < rvl) {last++; SWAP2ri(v[last],v[i],V[last],V[i],rtmp,itmp);}
12039f41f7fSStefano Zampini   }
12139f41f7fSStefano Zampini   SWAP2ri(v[0],v[last],V[0],V[last],rtmp,itmp);
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortRealWithArrayInt_Private(v,V,last-1));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortRealWithArrayInt_Private(v+last+1,V+last+1,right-(last+1)));
12439f41f7fSStefano Zampini   PetscFunctionReturn(0);
12539f41f7fSStefano Zampini }
12639f41f7fSStefano Zampini /*@
12739f41f7fSStefano Zampini    PetscSortRealWithArrayInt - Sorts an array of PetscReal in place in increasing order;
12839f41f7fSStefano Zampini        changes a second PetscInt array to match the sorted first array.
12939f41f7fSStefano Zampini 
13039f41f7fSStefano Zampini    Not Collective
13139f41f7fSStefano Zampini 
13239f41f7fSStefano Zampini    Input Parameters:
13339f41f7fSStefano Zampini +  n  - number of values
13439f41f7fSStefano Zampini .  i  - array of integers
13539f41f7fSStefano Zampini -  I - second array of integers
13639f41f7fSStefano Zampini 
13739f41f7fSStefano Zampini    Level: intermediate
13839f41f7fSStefano Zampini 
13939f41f7fSStefano Zampini .seealso: PetscSortReal()
14039f41f7fSStefano Zampini @*/
14139f41f7fSStefano Zampini PetscErrorCode  PetscSortRealWithArrayInt(PetscInt n,PetscReal r[],PetscInt Ii[])
14239f41f7fSStefano Zampini {
14339f41f7fSStefano Zampini   PetscInt       j,k,itmp;
14439f41f7fSStefano Zampini   PetscReal      rk,rtmp;
14539f41f7fSStefano Zampini 
14639f41f7fSStefano Zampini   PetscFunctionBegin;
147*dadcf809SJacob Faibussowitsch   PetscValidRealPointer(r,2);
148*dadcf809SJacob Faibussowitsch   PetscValidIntPointer(Ii,3);
14939f41f7fSStefano Zampini   if (n<8) {
15039f41f7fSStefano Zampini     for (k=0; k<n; k++) {
15139f41f7fSStefano Zampini       rk = r[k];
15239f41f7fSStefano Zampini       for (j=k+1; j<n; j++) {
15339f41f7fSStefano Zampini         if (rk > r[j]) {
15439f41f7fSStefano Zampini           SWAP2ri(r[k],r[j],Ii[k],Ii[j],rtmp,itmp);
15539f41f7fSStefano Zampini           rk = r[k];
15639f41f7fSStefano Zampini         }
15739f41f7fSStefano Zampini       }
15839f41f7fSStefano Zampini     }
15939f41f7fSStefano Zampini   } else {
1605f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSortRealWithArrayInt_Private(r,Ii,n-1));
16139f41f7fSStefano Zampini   }
16239f41f7fSStefano Zampini   PetscFunctionReturn(0);
16339f41f7fSStefano Zampini }
16439f41f7fSStefano Zampini 
16539f41f7fSStefano Zampini /*@
16639f41f7fSStefano Zampini   PetscFindReal - Finds a PetscReal in a sorted array of PetscReals
16739f41f7fSStefano Zampini 
16839f41f7fSStefano Zampini    Not Collective
16939f41f7fSStefano Zampini 
17039f41f7fSStefano Zampini    Input Parameters:
17139f41f7fSStefano Zampini +  key - the value to locate
17239f41f7fSStefano Zampini .  n   - number of values in the array
17339f41f7fSStefano Zampini .  ii  - array of values
17439f41f7fSStefano Zampini -  eps - tolerance used to compare
17539f41f7fSStefano Zampini 
17639f41f7fSStefano Zampini    Output Parameter:
17739f41f7fSStefano Zampini .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
17839f41f7fSStefano Zampini 
17939f41f7fSStefano Zampini    Level: intermediate
18039f41f7fSStefano Zampini 
18139f41f7fSStefano Zampini .seealso: PetscSortReal(), PetscSortRealWithArrayInt()
18239f41f7fSStefano Zampini @*/
18339f41f7fSStefano Zampini PetscErrorCode PetscFindReal(PetscReal key, PetscInt n, const PetscReal t[], PetscReal eps, PetscInt *loc)
18439f41f7fSStefano Zampini {
18539f41f7fSStefano Zampini   PetscInt lo = 0,hi = n;
18639f41f7fSStefano Zampini 
18739f41f7fSStefano Zampini   PetscFunctionBegin;
188*dadcf809SJacob Faibussowitsch   PetscValidIntPointer(loc,5);
18939f41f7fSStefano Zampini   if (!n) {*loc = -1; PetscFunctionReturn(0);}
190*dadcf809SJacob Faibussowitsch   PetscValidRealPointer(t,3);
1916a7c706bSVaclav Hapla   PetscCheckSorted(n,t);
19239f41f7fSStefano Zampini   while (hi - lo > 1) {
19339f41f7fSStefano Zampini     PetscInt mid = lo + (hi - lo)/2;
19439f41f7fSStefano Zampini     if (key < t[mid]) hi = mid;
19539f41f7fSStefano Zampini     else              lo = mid;
19639f41f7fSStefano Zampini   }
19739f41f7fSStefano Zampini   *loc = (PetscAbsReal(key - t[lo]) < eps) ? lo : -(lo + (key > t[lo]) + 1);
19839f41f7fSStefano Zampini   PetscFunctionReturn(0);
19939f41f7fSStefano Zampini }
200745b41b2SMatthew G. Knepley 
201745b41b2SMatthew G. Knepley /*@
202745b41b2SMatthew G. Knepley    PetscSortRemoveDupsReal - Sorts an array of doubles in place in increasing order removes all duplicate entries
203745b41b2SMatthew G. Knepley 
204745b41b2SMatthew G. Knepley    Not Collective
205745b41b2SMatthew G. Knepley 
206745b41b2SMatthew G. Knepley    Input Parameters:
207745b41b2SMatthew G. Knepley +  n  - number of values
208745b41b2SMatthew G. Knepley -  v  - array of doubles
209745b41b2SMatthew G. Knepley 
210745b41b2SMatthew G. Knepley    Output Parameter:
211745b41b2SMatthew G. Knepley .  n - number of non-redundant values
212745b41b2SMatthew G. Knepley 
213745b41b2SMatthew G. Knepley    Level: intermediate
214745b41b2SMatthew G. Knepley 
215745b41b2SMatthew G. Knepley .seealso: PetscSortReal(), PetscSortRemoveDupsInt()
216745b41b2SMatthew G. Knepley @*/
217745b41b2SMatthew G. Knepley PetscErrorCode  PetscSortRemoveDupsReal(PetscInt *n,PetscReal v[])
218745b41b2SMatthew G. Knepley {
219745b41b2SMatthew G. Knepley   PetscInt       i,s = 0,N = *n, b = 0;
220745b41b2SMatthew G. Knepley 
221745b41b2SMatthew G. Knepley   PetscFunctionBegin;
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortReal(N,v));
223745b41b2SMatthew G. Knepley   for (i=0; i<N-1; i++) {
224745b41b2SMatthew G. Knepley     if (v[b+s+1] != v[b]) {
225745b41b2SMatthew G. Knepley       v[b+1] = v[b+s+1]; b++;
226745b41b2SMatthew G. Knepley     } else s++;
227745b41b2SMatthew G. Knepley   }
228745b41b2SMatthew G. Knepley   *n = N - s;
229745b41b2SMatthew G. Knepley   PetscFunctionReturn(0);
230745b41b2SMatthew G. Knepley }
231745b41b2SMatthew G. Knepley 
232d97c5584SHong Zhang /*@
233021822bcSHong Zhang    PetscSortSplit - Quick-sort split of an array of PetscScalars in place.
234d97c5584SHong Zhang 
235d97c5584SHong Zhang    Not Collective
236d97c5584SHong Zhang 
237d97c5584SHong Zhang    Input Parameters:
238d97c5584SHong Zhang +  ncut  - splitig index
2396b867d5aSJose E. Roman -  n     - number of values to sort
240d97c5584SHong Zhang 
2416b867d5aSJose E. Roman    Input/Output Parameters:
2426b867d5aSJose E. Roman +  a     - array of values, on output the values are permuted such that its elements satisfy:
243d97c5584SHong Zhang            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
244d97c5584SHong Zhang            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
2456b867d5aSJose E. Roman -  idx   - index for array a, on output permuted accordingly
246d97c5584SHong Zhang 
247d97c5584SHong Zhang    Level: intermediate
248d97c5584SHong Zhang 
249d97c5584SHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation()
250d97c5584SHong Zhang @*/
2517087cfbeSBarry Smith PetscErrorCode  PetscSortSplit(PetscInt ncut,PetscInt n,PetscScalar a[],PetscInt idx[])
252d97c5584SHong Zhang {
253d97c5584SHong Zhang   PetscInt    i,mid,last,itmp,j,first;
254d97c5584SHong Zhang   PetscScalar d,tmp;
255d97c5584SHong Zhang   PetscReal   abskey;
256d97c5584SHong Zhang 
257d97c5584SHong Zhang   PetscFunctionBegin;
258d97c5584SHong Zhang   first = 0;
259d97c5584SHong Zhang   last  = n-1;
260d97c5584SHong Zhang   if (ncut < first || ncut > last) PetscFunctionReturn(0);
261d97c5584SHong Zhang 
262d97c5584SHong Zhang   while (1) {
263d97c5584SHong Zhang     mid    = first;
2642a274a78SSatish Balay     d      = a[mid];
2652a274a78SSatish Balay     abskey = PetscAbsScalar(d);
266d97c5584SHong Zhang     i      = last;
267d97c5584SHong Zhang     for (j = first + 1; j <= i; ++j) {
2682a274a78SSatish Balay       d = a[j];
2692a274a78SSatish Balay       if (PetscAbsScalar(d) >= abskey) {
270d97c5584SHong Zhang         ++mid;
271d97c5584SHong Zhang         /* interchange */
272d97c5584SHong Zhang         tmp = a[mid];  itmp = idx[mid];
273d97c5584SHong Zhang         a[mid] = a[j]; idx[mid] = idx[j];
274d97c5584SHong Zhang         a[j] = tmp;    idx[j] = itmp;
275d97c5584SHong Zhang       }
276d97c5584SHong Zhang     }
277d97c5584SHong Zhang 
278d97c5584SHong Zhang     /* interchange */
279d97c5584SHong Zhang     tmp = a[mid];      itmp = idx[mid];
280d97c5584SHong Zhang     a[mid] = a[first]; idx[mid] = idx[first];
281d97c5584SHong Zhang     a[first] = tmp;    idx[first] = itmp;
282d97c5584SHong Zhang 
283d97c5584SHong Zhang     /* test for while loop */
284a297a907SKarl Rupp     if (mid == ncut) break;
285a297a907SKarl Rupp     else if (mid > ncut) last = mid - 1;
286a297a907SKarl Rupp     else first = mid + 1;
287d97c5584SHong Zhang   }
288d97c5584SHong Zhang   PetscFunctionReturn(0);
289d97c5584SHong Zhang }
290021822bcSHong Zhang 
291021822bcSHong Zhang /*@
292021822bcSHong Zhang    PetscSortSplitReal - Quick-sort split of an array of PetscReals in place.
293021822bcSHong Zhang 
294021822bcSHong Zhang    Not Collective
295021822bcSHong Zhang 
296021822bcSHong Zhang    Input Parameters:
297021822bcSHong Zhang +  ncut  - splitig index
2986b867d5aSJose E. Roman -  n     - number of values to sort
299021822bcSHong Zhang 
3006b867d5aSJose E. Roman    Input/Output Parameters:
3016b867d5aSJose E. Roman +  a     - array of values, on output the values are permuted such that its elements satisfy:
302021822bcSHong Zhang            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
303021822bcSHong Zhang            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
3046b867d5aSJose E. Roman -  idx   - index for array a, on output permuted accordingly
305021822bcSHong Zhang 
306021822bcSHong Zhang    Level: intermediate
307021822bcSHong Zhang 
308021822bcSHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation()
309021822bcSHong Zhang @*/
3107087cfbeSBarry Smith PetscErrorCode  PetscSortSplitReal(PetscInt ncut,PetscInt n,PetscReal a[],PetscInt idx[])
311021822bcSHong Zhang {
312021822bcSHong Zhang   PetscInt  i,mid,last,itmp,j,first;
313021822bcSHong Zhang   PetscReal d,tmp;
314021822bcSHong Zhang   PetscReal abskey;
315021822bcSHong Zhang 
316021822bcSHong Zhang   PetscFunctionBegin;
317021822bcSHong Zhang   first = 0;
318021822bcSHong Zhang   last  = n-1;
319021822bcSHong Zhang   if (ncut < first || ncut > last) PetscFunctionReturn(0);
320021822bcSHong Zhang 
321021822bcSHong Zhang   while (1) {
322021822bcSHong Zhang     mid    = first;
3232a274a78SSatish Balay     d      = a[mid];
3242a274a78SSatish Balay     abskey = PetscAbsReal(d);
325021822bcSHong Zhang     i      = last;
326021822bcSHong Zhang     for (j = first + 1; j <= i; ++j) {
3272a274a78SSatish Balay       d = a[j];
3282a274a78SSatish Balay       if (PetscAbsReal(d) >= abskey) {
329021822bcSHong Zhang         ++mid;
330021822bcSHong Zhang         /* interchange */
331021822bcSHong Zhang         tmp = a[mid];  itmp = idx[mid];
332021822bcSHong Zhang         a[mid] = a[j]; idx[mid] = idx[j];
333021822bcSHong Zhang         a[j] = tmp;    idx[j] = itmp;
334021822bcSHong Zhang       }
335021822bcSHong Zhang     }
336021822bcSHong Zhang 
337021822bcSHong Zhang     /* interchange */
338021822bcSHong Zhang     tmp = a[mid];      itmp = idx[mid];
339021822bcSHong Zhang     a[mid] = a[first]; idx[mid] = idx[first];
340021822bcSHong Zhang     a[first] = tmp;    idx[first] = itmp;
341021822bcSHong Zhang 
342021822bcSHong Zhang     /* test for while loop */
343a297a907SKarl Rupp     if (mid == ncut) break;
344a297a907SKarl Rupp     else if (mid > ncut) last = mid - 1;
345a297a907SKarl Rupp     else first = mid + 1;
346021822bcSHong Zhang   }
347021822bcSHong Zhang   PetscFunctionReturn(0);
348021822bcSHong Zhang }
349