xref: /petsc/src/sys/utils/sortd.c (revision 39f41f7f10bf62dde9956fb6ea07e97418fcde10)
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*/
9*39f41f7fSStefano Zampini #include <petsc/private/petscimpl.h>
10e5c89e4eSSatish Balay 
11e5c89e4eSSatish Balay #define SWAP(a,b,t) {t=a;a=b;b=t;}
12e5c89e4eSSatish Balay 
13e5c89e4eSSatish Balay /* A simple version of quicksort; taken from Kernighan and Ritchie, page 87 */
14e5c89e4eSSatish Balay static PetscErrorCode PetscSortReal_Private(PetscReal *v,PetscInt right)
15e5c89e4eSSatish Balay {
16e5c89e4eSSatish Balay   PetscInt  i,last;
17e5c89e4eSSatish Balay   PetscReal vl,tmp;
18e5c89e4eSSatish Balay 
19e5c89e4eSSatish Balay   PetscFunctionBegin;
20e5c89e4eSSatish Balay   if (right <= 1) {
21e5c89e4eSSatish Balay     if (right == 1) {
22e5c89e4eSSatish Balay       if (v[0] > v[1]) SWAP(v[0],v[1],tmp);
23e5c89e4eSSatish Balay     }
24e5c89e4eSSatish Balay     PetscFunctionReturn(0);
25e5c89e4eSSatish Balay   }
26e5c89e4eSSatish Balay   SWAP(v[0],v[right/2],tmp);
27e5c89e4eSSatish Balay   vl   = v[0];
28e5c89e4eSSatish Balay   last = 0;
29e5c89e4eSSatish Balay   for (i=1; i<=right; i++) {
30e5c89e4eSSatish Balay     if (v[i] < vl) {last++; SWAP(v[last],v[i],tmp);}
31e5c89e4eSSatish Balay   }
32e5c89e4eSSatish Balay   SWAP(v[0],v[last],tmp);
33e5c89e4eSSatish Balay   PetscSortReal_Private(v,last-1);
34e5c89e4eSSatish Balay   PetscSortReal_Private(v+last+1,right-(last+1));
35e5c89e4eSSatish Balay   PetscFunctionReturn(0);
36e5c89e4eSSatish Balay }
37e5c89e4eSSatish Balay 
38e5c89e4eSSatish Balay /*@
39e5c89e4eSSatish Balay    PetscSortReal - Sorts an array of doubles in place in increasing order.
40e5c89e4eSSatish Balay 
41e5c89e4eSSatish Balay    Not Collective
42e5c89e4eSSatish Balay 
43e5c89e4eSSatish Balay    Input Parameters:
44e5c89e4eSSatish Balay +  n  - number of values
45e5c89e4eSSatish Balay -  v  - array of doubles
46e5c89e4eSSatish Balay 
47e5c89e4eSSatish Balay    Level: intermediate
48e5c89e4eSSatish Balay 
49e5c89e4eSSatish Balay    Concepts: sorting^doubles
50e5c89e4eSSatish Balay 
51*39f41f7fSStefano Zampini .seealso: PetscSortInt(), PetscSortRealWithPermutation(), PetscSortRealWithArrayInt()
52e5c89e4eSSatish Balay @*/
537087cfbeSBarry Smith PetscErrorCode  PetscSortReal(PetscInt n,PetscReal v[])
54e5c89e4eSSatish Balay {
55e5c89e4eSSatish Balay   PetscInt  j,k;
56e5c89e4eSSatish Balay   PetscReal tmp,vk;
57e5c89e4eSSatish Balay 
58e5c89e4eSSatish Balay   PetscFunctionBegin;
59*39f41f7fSStefano Zampini   PetscValidPointer(v,2);
60e5c89e4eSSatish Balay   if (n<8) {
61e5c89e4eSSatish Balay     for (k=0; k<n; k++) {
62e5c89e4eSSatish Balay       vk = v[k];
63e5c89e4eSSatish Balay       for (j=k+1; j<n; j++) {
64e5c89e4eSSatish Balay         if (vk > v[j]) {
65e5c89e4eSSatish Balay           SWAP(v[k],v[j],tmp);
66e5c89e4eSSatish Balay           vk = v[k];
67e5c89e4eSSatish Balay         }
68e5c89e4eSSatish Balay       }
69e5c89e4eSSatish Balay     }
70a297a907SKarl Rupp   } else PetscSortReal_Private(v,n-1);
71e5c89e4eSSatish Balay   PetscFunctionReturn(0);
72e5c89e4eSSatish Balay }
73e5c89e4eSSatish Balay 
74*39f41f7fSStefano Zampini #define SWAP2ri(a,b,c,d,rt,it) {rt=a;a=b;b=rt;it=c;c=d;d=it;}
75*39f41f7fSStefano Zampini 
76*39f41f7fSStefano Zampini /* modified from PetscSortIntWithArray_Private */
77*39f41f7fSStefano Zampini static PetscErrorCode PetscSortRealWithArrayInt_Private(PetscReal *v,PetscInt *V,PetscInt right)
78*39f41f7fSStefano Zampini {
79*39f41f7fSStefano Zampini   PetscErrorCode ierr;
80*39f41f7fSStefano Zampini   PetscInt       i,last,itmp;
81*39f41f7fSStefano Zampini   PetscReal      rvl,rtmp;
82*39f41f7fSStefano Zampini 
83*39f41f7fSStefano Zampini   PetscFunctionBegin;
84*39f41f7fSStefano Zampini   if (right <= 1) {
85*39f41f7fSStefano Zampini     if (right == 1) {
86*39f41f7fSStefano Zampini       if (v[0] > v[1]) SWAP2ri(v[0],v[1],V[0],V[1],rtmp,itmp);
87*39f41f7fSStefano Zampini     }
88*39f41f7fSStefano Zampini     PetscFunctionReturn(0);
89*39f41f7fSStefano Zampini   }
90*39f41f7fSStefano Zampini   SWAP2ri(v[0],v[right/2],V[0],V[right/2],rtmp,itmp);
91*39f41f7fSStefano Zampini   rvl  = v[0];
92*39f41f7fSStefano Zampini   last = 0;
93*39f41f7fSStefano Zampini   for (i=1; i<=right; i++) {
94*39f41f7fSStefano Zampini     if (v[i] < rvl) {last++; SWAP2ri(v[last],v[i],V[last],V[i],rtmp,itmp);}
95*39f41f7fSStefano Zampini   }
96*39f41f7fSStefano Zampini   SWAP2ri(v[0],v[last],V[0],V[last],rtmp,itmp);
97*39f41f7fSStefano Zampini   ierr = PetscSortRealWithArrayInt_Private(v,V,last-1);CHKERRQ(ierr);
98*39f41f7fSStefano Zampini   ierr = PetscSortRealWithArrayInt_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr);
99*39f41f7fSStefano Zampini   PetscFunctionReturn(0);
100*39f41f7fSStefano Zampini }
101*39f41f7fSStefano Zampini /*@
102*39f41f7fSStefano Zampini    PetscSortRealWithArrayInt - Sorts an array of PetscReal in place in increasing order;
103*39f41f7fSStefano Zampini        changes a second PetscInt array to match the sorted first array.
104*39f41f7fSStefano Zampini 
105*39f41f7fSStefano Zampini    Not Collective
106*39f41f7fSStefano Zampini 
107*39f41f7fSStefano Zampini    Input Parameters:
108*39f41f7fSStefano Zampini +  n  - number of values
109*39f41f7fSStefano Zampini .  i  - array of integers
110*39f41f7fSStefano Zampini -  I - second array of integers
111*39f41f7fSStefano Zampini 
112*39f41f7fSStefano Zampini    Level: intermediate
113*39f41f7fSStefano Zampini 
114*39f41f7fSStefano Zampini    Concepts: sorting^ints with array
115*39f41f7fSStefano Zampini 
116*39f41f7fSStefano Zampini .seealso: PetscSortReal()
117*39f41f7fSStefano Zampini @*/
118*39f41f7fSStefano Zampini PetscErrorCode  PetscSortRealWithArrayInt(PetscInt n,PetscReal r[],PetscInt Ii[])
119*39f41f7fSStefano Zampini {
120*39f41f7fSStefano Zampini   PetscErrorCode ierr;
121*39f41f7fSStefano Zampini   PetscInt       j,k,itmp;
122*39f41f7fSStefano Zampini   PetscReal      rk,rtmp;
123*39f41f7fSStefano Zampini 
124*39f41f7fSStefano Zampini   PetscFunctionBegin;
125*39f41f7fSStefano Zampini   PetscValidPointer(r,2);
126*39f41f7fSStefano Zampini   PetscValidPointer(Ii,3);
127*39f41f7fSStefano Zampini   if (n<8) {
128*39f41f7fSStefano Zampini     for (k=0; k<n; k++) {
129*39f41f7fSStefano Zampini       rk = r[k];
130*39f41f7fSStefano Zampini       for (j=k+1; j<n; j++) {
131*39f41f7fSStefano Zampini         if (rk > r[j]) {
132*39f41f7fSStefano Zampini           SWAP2ri(r[k],r[j],Ii[k],Ii[j],rtmp,itmp);
133*39f41f7fSStefano Zampini           rk = r[k];
134*39f41f7fSStefano Zampini         }
135*39f41f7fSStefano Zampini       }
136*39f41f7fSStefano Zampini     }
137*39f41f7fSStefano Zampini   } else {
138*39f41f7fSStefano Zampini     ierr = PetscSortRealWithArrayInt_Private(r,Ii,n-1);CHKERRQ(ierr);
139*39f41f7fSStefano Zampini   }
140*39f41f7fSStefano Zampini   PetscFunctionReturn(0);
141*39f41f7fSStefano Zampini }
142*39f41f7fSStefano Zampini 
143*39f41f7fSStefano Zampini /*@
144*39f41f7fSStefano Zampini   PetscFindReal - Finds a PetscReal in a sorted array of PetscReals
145*39f41f7fSStefano Zampini 
146*39f41f7fSStefano Zampini    Not Collective
147*39f41f7fSStefano Zampini 
148*39f41f7fSStefano Zampini    Input Parameters:
149*39f41f7fSStefano Zampini +  key - the value to locate
150*39f41f7fSStefano Zampini .  n   - number of values in the array
151*39f41f7fSStefano Zampini .  ii  - array of values
152*39f41f7fSStefano Zampini -  eps - tolerance used to compare
153*39f41f7fSStefano Zampini 
154*39f41f7fSStefano Zampini    Output Parameter:
155*39f41f7fSStefano Zampini .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
156*39f41f7fSStefano Zampini 
157*39f41f7fSStefano Zampini    Level: intermediate
158*39f41f7fSStefano Zampini 
159*39f41f7fSStefano Zampini    Concepts: sorting^ints
160*39f41f7fSStefano Zampini 
161*39f41f7fSStefano Zampini .seealso: PetscSortReal(), PetscSortRealWithArrayInt()
162*39f41f7fSStefano Zampini @*/
163*39f41f7fSStefano Zampini PetscErrorCode PetscFindReal(PetscReal key, PetscInt n, const PetscReal t[], PetscReal eps, PetscInt *loc)
164*39f41f7fSStefano Zampini {
165*39f41f7fSStefano Zampini   PetscInt lo = 0,hi = n;
166*39f41f7fSStefano Zampini 
167*39f41f7fSStefano Zampini   PetscFunctionBegin;
168*39f41f7fSStefano Zampini   PetscValidPointer(loc,4);
169*39f41f7fSStefano Zampini   if (!n) {*loc = -1; PetscFunctionReturn(0);}
170*39f41f7fSStefano Zampini   PetscValidPointer(t,3);
171*39f41f7fSStefano Zampini   while (hi - lo > 1) {
172*39f41f7fSStefano Zampini     PetscInt mid = lo + (hi - lo)/2;
173*39f41f7fSStefano Zampini     if (key < t[mid]) hi = mid;
174*39f41f7fSStefano Zampini     else              lo = mid;
175*39f41f7fSStefano Zampini   }
176*39f41f7fSStefano Zampini   *loc = (PetscAbsReal(key - t[lo]) < eps) ? lo : -(lo + (key > t[lo]) + 1);
177*39f41f7fSStefano Zampini   PetscFunctionReturn(0);
178*39f41f7fSStefano Zampini }
179745b41b2SMatthew G. Knepley 
180745b41b2SMatthew G. Knepley /*@
181745b41b2SMatthew G. Knepley    PetscSortRemoveDupsReal - Sorts an array of doubles in place in increasing order removes all duplicate entries
182745b41b2SMatthew G. Knepley 
183745b41b2SMatthew G. Knepley    Not Collective
184745b41b2SMatthew G. Knepley 
185745b41b2SMatthew G. Knepley    Input Parameters:
186745b41b2SMatthew G. Knepley +  n  - number of values
187745b41b2SMatthew G. Knepley -  v  - array of doubles
188745b41b2SMatthew G. Knepley 
189745b41b2SMatthew G. Knepley    Output Parameter:
190745b41b2SMatthew G. Knepley .  n - number of non-redundant values
191745b41b2SMatthew G. Knepley 
192745b41b2SMatthew G. Knepley    Level: intermediate
193745b41b2SMatthew G. Knepley 
194745b41b2SMatthew G. Knepley    Concepts: sorting^doubles
195745b41b2SMatthew G. Knepley 
196745b41b2SMatthew G. Knepley .seealso: PetscSortReal(), PetscSortRemoveDupsInt()
197745b41b2SMatthew G. Knepley @*/
198745b41b2SMatthew G. Knepley PetscErrorCode  PetscSortRemoveDupsReal(PetscInt *n,PetscReal v[])
199745b41b2SMatthew G. Knepley {
200745b41b2SMatthew G. Knepley   PetscErrorCode ierr;
201745b41b2SMatthew G. Knepley   PetscInt       i,s = 0,N = *n, b = 0;
202745b41b2SMatthew G. Knepley 
203745b41b2SMatthew G. Knepley   PetscFunctionBegin;
204745b41b2SMatthew G. Knepley   ierr = PetscSortReal(N,v);CHKERRQ(ierr);
205745b41b2SMatthew G. Knepley   for (i=0; i<N-1; i++) {
206745b41b2SMatthew G. Knepley     if (v[b+s+1] != v[b]) {
207745b41b2SMatthew G. Knepley       v[b+1] = v[b+s+1]; b++;
208745b41b2SMatthew G. Knepley     } else s++;
209745b41b2SMatthew G. Knepley   }
210745b41b2SMatthew G. Knepley   *n = N - s;
211745b41b2SMatthew G. Knepley   PetscFunctionReturn(0);
212745b41b2SMatthew G. Knepley }
213745b41b2SMatthew G. Knepley 
214d97c5584SHong Zhang /*@
215021822bcSHong Zhang    PetscSortSplit - Quick-sort split of an array of PetscScalars in place.
216d97c5584SHong Zhang 
217d97c5584SHong Zhang    Not Collective
218d97c5584SHong Zhang 
219d97c5584SHong Zhang    Input Parameters:
220d97c5584SHong Zhang +  ncut  - splitig index
221d97c5584SHong Zhang .  n     - number of values to sort
222d97c5584SHong Zhang .  a     - array of values
223d97c5584SHong Zhang -  idx   - index for array a
224d97c5584SHong Zhang 
225d97c5584SHong Zhang    Output Parameters:
226d97c5584SHong Zhang +  a     - permuted array of values such that its elements satisfy:
227d97c5584SHong Zhang            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
228d97c5584SHong Zhang            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
229d97c5584SHong Zhang -  idx   - permuted index of array a
230d97c5584SHong Zhang 
231d97c5584SHong Zhang    Level: intermediate
232d97c5584SHong Zhang 
233d97c5584SHong Zhang    Concepts: sorting^doubles
234d97c5584SHong Zhang 
235d97c5584SHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation()
236d97c5584SHong Zhang @*/
2377087cfbeSBarry Smith PetscErrorCode  PetscSortSplit(PetscInt ncut,PetscInt n,PetscScalar a[],PetscInt idx[])
238d97c5584SHong Zhang {
239d97c5584SHong Zhang   PetscInt    i,mid,last,itmp,j,first;
240d97c5584SHong Zhang   PetscScalar d,tmp;
241d97c5584SHong Zhang   PetscReal   abskey;
242d97c5584SHong Zhang 
243d97c5584SHong Zhang   PetscFunctionBegin;
244d97c5584SHong Zhang   first = 0;
245d97c5584SHong Zhang   last  = n-1;
246d97c5584SHong Zhang   if (ncut < first || ncut > last) PetscFunctionReturn(0);
247d97c5584SHong Zhang 
248d97c5584SHong Zhang   while (1) {
249d97c5584SHong Zhang     mid    = first;
2502a274a78SSatish Balay     d      = a[mid];
2512a274a78SSatish Balay     abskey = PetscAbsScalar(d);
252d97c5584SHong Zhang     i      = last;
253d97c5584SHong Zhang     for (j = first + 1; j <= i; ++j) {
2542a274a78SSatish Balay       d = a[j];
2552a274a78SSatish Balay       if (PetscAbsScalar(d) >= abskey) {
256d97c5584SHong Zhang         ++mid;
257d97c5584SHong Zhang         /* interchange */
258d97c5584SHong Zhang         tmp = a[mid];  itmp = idx[mid];
259d97c5584SHong Zhang         a[mid] = a[j]; idx[mid] = idx[j];
260d97c5584SHong Zhang         a[j] = tmp;    idx[j] = itmp;
261d97c5584SHong Zhang       }
262d97c5584SHong Zhang     }
263d97c5584SHong Zhang 
264d97c5584SHong Zhang     /* interchange */
265d97c5584SHong Zhang     tmp = a[mid];      itmp = idx[mid];
266d97c5584SHong Zhang     a[mid] = a[first]; idx[mid] = idx[first];
267d97c5584SHong Zhang     a[first] = tmp;    idx[first] = itmp;
268d97c5584SHong Zhang 
269d97c5584SHong Zhang     /* test for while loop */
270a297a907SKarl Rupp     if (mid == ncut) break;
271a297a907SKarl Rupp     else if (mid > ncut) last = mid - 1;
272a297a907SKarl Rupp     else first = mid + 1;
273d97c5584SHong Zhang   }
274d97c5584SHong Zhang   PetscFunctionReturn(0);
275d97c5584SHong Zhang }
276021822bcSHong Zhang 
277021822bcSHong Zhang /*@
278021822bcSHong Zhang    PetscSortSplitReal - Quick-sort split of an array of PetscReals in place.
279021822bcSHong Zhang 
280021822bcSHong Zhang    Not Collective
281021822bcSHong Zhang 
282021822bcSHong Zhang    Input Parameters:
283021822bcSHong Zhang +  ncut  - splitig index
284021822bcSHong Zhang .  n     - number of values to sort
285021822bcSHong Zhang .  a     - array of values in PetscReal
286021822bcSHong Zhang -  idx   - index for array a
287021822bcSHong Zhang 
288021822bcSHong Zhang    Output Parameters:
289021822bcSHong Zhang +  a     - permuted array of real values such that its elements satisfy:
290021822bcSHong Zhang            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
291021822bcSHong Zhang            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
292021822bcSHong Zhang -  idx   - permuted index of array a
293021822bcSHong Zhang 
294021822bcSHong Zhang    Level: intermediate
295021822bcSHong Zhang 
296021822bcSHong Zhang    Concepts: sorting^doubles
297021822bcSHong Zhang 
298021822bcSHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation()
299021822bcSHong Zhang @*/
3007087cfbeSBarry Smith PetscErrorCode  PetscSortSplitReal(PetscInt ncut,PetscInt n,PetscReal a[],PetscInt idx[])
301021822bcSHong Zhang {
302021822bcSHong Zhang   PetscInt  i,mid,last,itmp,j,first;
303021822bcSHong Zhang   PetscReal d,tmp;
304021822bcSHong Zhang   PetscReal abskey;
305021822bcSHong Zhang 
306021822bcSHong Zhang   PetscFunctionBegin;
307021822bcSHong Zhang   first = 0;
308021822bcSHong Zhang   last  = n-1;
309021822bcSHong Zhang   if (ncut < first || ncut > last) PetscFunctionReturn(0);
310021822bcSHong Zhang 
311021822bcSHong Zhang   while (1) {
312021822bcSHong Zhang     mid    = first;
3132a274a78SSatish Balay     d      = a[mid];
3142a274a78SSatish Balay     abskey = PetscAbsReal(d);
315021822bcSHong Zhang     i      = last;
316021822bcSHong Zhang     for (j = first + 1; j <= i; ++j) {
3172a274a78SSatish Balay       d = a[j];
3182a274a78SSatish Balay       if (PetscAbsReal(d) >= abskey) {
319021822bcSHong Zhang         ++mid;
320021822bcSHong Zhang         /* interchange */
321021822bcSHong Zhang         tmp = a[mid];  itmp = idx[mid];
322021822bcSHong Zhang         a[mid] = a[j]; idx[mid] = idx[j];
323021822bcSHong Zhang         a[j] = tmp;    idx[j] = itmp;
324021822bcSHong Zhang       }
325021822bcSHong Zhang     }
326021822bcSHong Zhang 
327021822bcSHong Zhang     /* interchange */
328021822bcSHong Zhang     tmp = a[mid];      itmp = idx[mid];
329021822bcSHong Zhang     a[mid] = a[first]; idx[mid] = idx[first];
330021822bcSHong Zhang     a[first] = tmp;    idx[first] = itmp;
331021822bcSHong Zhang 
332021822bcSHong Zhang     /* test for while loop */
333a297a907SKarl Rupp     if (mid == ncut) break;
334a297a907SKarl Rupp     else if (mid > ncut) last = mid - 1;
335a297a907SKarl Rupp     else first = mid + 1;
336021822bcSHong Zhang   }
337021822bcSHong Zhang   PetscFunctionReturn(0);
338021822bcSHong Zhang }
339021822bcSHong Zhang 
340