xref: /petsc/src/sys/utils/sortd.c (revision 7d0a6c19129e7069c8a40e210b34ed62989173db)
1*7d0a6c19SBarry 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  */
8d382aafbSBarry Smith #include "petscsys.h"                /*I  "petscsys.h"  I*/
9e5c89e4eSSatish Balay 
10e5c89e4eSSatish Balay #define SWAP(a,b,t) {t=a;a=b;b=t;}
11e5c89e4eSSatish Balay 
12e5c89e4eSSatish Balay #undef __FUNCT__
13e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortReal_Private"
14e5c89e4eSSatish Balay /* A simple version of quicksort; taken from Kernighan and Ritchie, page 87 */
15e5c89e4eSSatish Balay static PetscErrorCode PetscSortReal_Private(PetscReal *v,PetscInt right)
16e5c89e4eSSatish Balay {
17e5c89e4eSSatish Balay   PetscInt  i,last;
18e5c89e4eSSatish Balay   PetscReal vl,tmp;
19e5c89e4eSSatish Balay 
20e5c89e4eSSatish Balay   PetscFunctionBegin;
21e5c89e4eSSatish Balay   if (right <= 1) {
22e5c89e4eSSatish Balay     if (right == 1) {
23e5c89e4eSSatish Balay       if (v[0] > v[1]) SWAP(v[0],v[1],tmp);
24e5c89e4eSSatish Balay     }
25e5c89e4eSSatish Balay     PetscFunctionReturn(0);
26e5c89e4eSSatish Balay   }
27e5c89e4eSSatish Balay   SWAP(v[0],v[right/2],tmp);
28e5c89e4eSSatish Balay   vl   = v[0];
29e5c89e4eSSatish Balay   last = 0;
30e5c89e4eSSatish Balay   for (i=1; i<=right; i++) {
31e5c89e4eSSatish Balay     if (v[i] < vl) {last++; SWAP(v[last],v[i],tmp);}
32e5c89e4eSSatish Balay   }
33e5c89e4eSSatish Balay   SWAP(v[0],v[last],tmp);
34e5c89e4eSSatish Balay   PetscSortReal_Private(v,last-1);
35e5c89e4eSSatish Balay   PetscSortReal_Private(v+last+1,right-(last+1));
36e5c89e4eSSatish Balay   PetscFunctionReturn(0);
37e5c89e4eSSatish Balay }
38e5c89e4eSSatish Balay 
39e5c89e4eSSatish Balay #undef __FUNCT__
40e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortReal"
41e5c89e4eSSatish Balay /*@
42e5c89e4eSSatish Balay    PetscSortReal - Sorts an array of doubles in place in increasing order.
43e5c89e4eSSatish Balay 
44e5c89e4eSSatish Balay    Not Collective
45e5c89e4eSSatish Balay 
46e5c89e4eSSatish Balay    Input Parameters:
47e5c89e4eSSatish Balay +  n  - number of values
48e5c89e4eSSatish Balay -  v  - array of doubles
49e5c89e4eSSatish Balay 
50e5c89e4eSSatish Balay    Level: intermediate
51e5c89e4eSSatish Balay 
52e5c89e4eSSatish Balay    Concepts: sorting^doubles
53e5c89e4eSSatish Balay 
54e5c89e4eSSatish Balay .seealso: PetscSortInt(), PetscSortRealWithPermutation()
55e5c89e4eSSatish Balay @*/
567087cfbeSBarry Smith PetscErrorCode  PetscSortReal(PetscInt n,PetscReal v[])
57e5c89e4eSSatish Balay {
58e5c89e4eSSatish Balay   PetscInt  j,k;
59e5c89e4eSSatish Balay   PetscReal tmp,vk;
60e5c89e4eSSatish Balay 
61e5c89e4eSSatish Balay   PetscFunctionBegin;
62e5c89e4eSSatish Balay   if (n<8) {
63e5c89e4eSSatish Balay     for (k=0; k<n; k++) {
64e5c89e4eSSatish Balay       vk = v[k];
65e5c89e4eSSatish Balay       for (j=k+1; j<n; j++) {
66e5c89e4eSSatish Balay 	if (vk > v[j]) {
67e5c89e4eSSatish Balay 	  SWAP(v[k],v[j],tmp);
68e5c89e4eSSatish Balay 	  vk = v[k];
69e5c89e4eSSatish Balay 	}
70e5c89e4eSSatish Balay       }
71e5c89e4eSSatish Balay     }
72e5c89e4eSSatish Balay   } else {
73e5c89e4eSSatish Balay     PetscSortReal_Private(v,n-1);
74e5c89e4eSSatish Balay   }
75e5c89e4eSSatish Balay   PetscFunctionReturn(0);
76e5c89e4eSSatish Balay }
77e5c89e4eSSatish Balay 
78d97c5584SHong Zhang #undef __FUNCT__
79d97c5584SHong Zhang #define __FUNCT__ "PetscSortSplit"
80d97c5584SHong Zhang /*@
81021822bcSHong Zhang    PetscSortSplit - Quick-sort split of an array of PetscScalars in place.
82d97c5584SHong Zhang 
83d97c5584SHong Zhang    Not Collective
84d97c5584SHong Zhang 
85d97c5584SHong Zhang    Input Parameters:
86d97c5584SHong Zhang +  ncut  - splitig index
87d97c5584SHong Zhang .  n     - number of values to sort
88d97c5584SHong Zhang .  a     - array of values
89d97c5584SHong Zhang -  idx   - index for array a
90d97c5584SHong Zhang 
91d97c5584SHong Zhang    Output Parameters:
92d97c5584SHong Zhang +  a     - permuted array of values such that its elements satisfy:
93d97c5584SHong Zhang            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
94d97c5584SHong Zhang            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
95d97c5584SHong Zhang -  idx   - permuted index of array a
96d97c5584SHong Zhang 
97d97c5584SHong Zhang    Level: intermediate
98d97c5584SHong Zhang 
99d97c5584SHong Zhang    Concepts: sorting^doubles
100d97c5584SHong Zhang 
101d97c5584SHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation()
102d97c5584SHong Zhang @*/
1037087cfbeSBarry Smith PetscErrorCode  PetscSortSplit(PetscInt ncut,PetscInt n,PetscScalar a[],PetscInt idx[])
104d97c5584SHong Zhang {
105d97c5584SHong Zhang   PetscInt    i,mid,last,itmp,j,first;
106d97c5584SHong Zhang   PetscScalar d,tmp;
107d97c5584SHong Zhang   PetscReal   abskey;
108d97c5584SHong Zhang 
109d97c5584SHong Zhang   PetscFunctionBegin;
110d97c5584SHong Zhang   first = 0;
111d97c5584SHong Zhang   last  = n-1;
112d97c5584SHong Zhang   if (ncut < first || ncut > last) PetscFunctionReturn(0);
113d97c5584SHong Zhang 
114d97c5584SHong Zhang   while (1){
115d97c5584SHong Zhang     mid = first;
116d97c5584SHong Zhang     abskey = (d = a[mid],PetscAbsScalar(d));
117d97c5584SHong Zhang     i = last;
118d97c5584SHong Zhang     for (j = first + 1; j <= i; ++j) {
119d97c5584SHong Zhang       if ((d = a[j],PetscAbsScalar(d)) >= abskey) {
120d97c5584SHong Zhang         ++mid;
121d97c5584SHong Zhang         /* interchange */
122d97c5584SHong Zhang         tmp = a[mid];  itmp = idx[mid];
123d97c5584SHong Zhang         a[mid] = a[j]; idx[mid] = idx[j];
124d97c5584SHong Zhang         a[j] = tmp;    idx[j] = itmp;
125d97c5584SHong Zhang       }
126d97c5584SHong Zhang     }
127d97c5584SHong Zhang 
128d97c5584SHong Zhang     /* interchange */
129d97c5584SHong Zhang     tmp = a[mid];      itmp = idx[mid];
130d97c5584SHong Zhang     a[mid] = a[first]; idx[mid] = idx[first];
131d97c5584SHong Zhang     a[first] = tmp;    idx[first] = itmp;
132d97c5584SHong Zhang 
133d97c5584SHong Zhang     /* test for while loop */
134d97c5584SHong Zhang     if (mid == ncut) {
135d97c5584SHong Zhang       break;
136d97c5584SHong Zhang     } else if (mid > ncut){
137d97c5584SHong Zhang       last = mid - 1;
138d97c5584SHong Zhang     } else {
139d97c5584SHong Zhang       first = mid + 1;
140d97c5584SHong Zhang     }
141d97c5584SHong Zhang   }
142d97c5584SHong Zhang   PetscFunctionReturn(0);
143d97c5584SHong Zhang }
144021822bcSHong Zhang 
145021822bcSHong Zhang #undef __FUNCT__
146021822bcSHong Zhang #define __FUNCT__ "PetscSortSplitReal"
147021822bcSHong Zhang /*@
148021822bcSHong Zhang    PetscSortSplitReal - Quick-sort split of an array of PetscReals in place.
149021822bcSHong Zhang 
150021822bcSHong Zhang    Not Collective
151021822bcSHong Zhang 
152021822bcSHong Zhang    Input Parameters:
153021822bcSHong Zhang +  ncut  - splitig index
154021822bcSHong Zhang .  n     - number of values to sort
155021822bcSHong Zhang .  a     - array of values in PetscReal
156021822bcSHong Zhang -  idx   - index for array a
157021822bcSHong Zhang 
158021822bcSHong Zhang    Output Parameters:
159021822bcSHong Zhang +  a     - permuted array of real values such that its elements satisfy:
160021822bcSHong Zhang            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
161021822bcSHong Zhang            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
162021822bcSHong Zhang -  idx   - permuted index of array a
163021822bcSHong Zhang 
164021822bcSHong Zhang    Level: intermediate
165021822bcSHong Zhang 
166021822bcSHong Zhang    Concepts: sorting^doubles
167021822bcSHong Zhang 
168021822bcSHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation()
169021822bcSHong Zhang @*/
1707087cfbeSBarry Smith PetscErrorCode  PetscSortSplitReal(PetscInt ncut,PetscInt n,PetscReal a[],PetscInt idx[])
171021822bcSHong Zhang {
172021822bcSHong Zhang   PetscInt    i,mid,last,itmp,j,first;
173021822bcSHong Zhang   PetscReal   d,tmp;
174021822bcSHong Zhang   PetscReal   abskey;
175021822bcSHong Zhang 
176021822bcSHong Zhang   PetscFunctionBegin;
177021822bcSHong Zhang   first = 0;
178021822bcSHong Zhang   last  = n-1;
179021822bcSHong Zhang   if (ncut < first || ncut > last) PetscFunctionReturn(0);
180021822bcSHong Zhang 
181021822bcSHong Zhang   while (1){
182021822bcSHong Zhang     mid = first;
183021822bcSHong Zhang     abskey = (d = a[mid],PetscAbsScalar(d));
184021822bcSHong Zhang     i = last;
185021822bcSHong Zhang     for (j = first + 1; j <= i; ++j) {
186021822bcSHong Zhang       if ((d = a[j],PetscAbsScalar(d)) >= abskey) {
187021822bcSHong Zhang         ++mid;
188021822bcSHong Zhang         /* interchange */
189021822bcSHong Zhang         tmp = a[mid];  itmp = idx[mid];
190021822bcSHong Zhang         a[mid] = a[j]; idx[mid] = idx[j];
191021822bcSHong Zhang         a[j] = tmp;    idx[j] = itmp;
192021822bcSHong Zhang       }
193021822bcSHong Zhang     }
194021822bcSHong Zhang 
195021822bcSHong Zhang     /* interchange */
196021822bcSHong Zhang     tmp = a[mid];      itmp = idx[mid];
197021822bcSHong Zhang     a[mid] = a[first]; idx[mid] = idx[first];
198021822bcSHong Zhang     a[first] = tmp;    idx[first] = itmp;
199021822bcSHong Zhang 
200021822bcSHong Zhang     /* test for while loop */
201021822bcSHong Zhang     if (mid == ncut) {
202021822bcSHong Zhang       break;
203021822bcSHong Zhang     } else if (mid > ncut){
204021822bcSHong Zhang       last = mid - 1;
205021822bcSHong Zhang     } else {
206021822bcSHong Zhang       first = mid + 1;
207021822bcSHong Zhang     }
208021822bcSHong Zhang   }
209021822bcSHong Zhang   PetscFunctionReturn(0);
210021822bcSHong Zhang }
211021822bcSHong Zhang 
212