Actual source code: sortd.c

petsc-master 2016-09-25
Report Typos and Errors
  2: /*
  3:    This file contains routines for sorting doubles.  Values are sorted in place.
  4:    These are provided because the general sort routines incur a great deal
  5:    of overhead in calling the comparision routines.

  7:  */
 8:  #include <petscsys.h>

 10: #define SWAP(a,b,t) {t=a;a=b;b=t;}

 14: /* A simple version of quicksort; taken from Kernighan and Ritchie, page 87 */
 15: static PetscErrorCode PetscSortReal_Private(PetscReal *v,PetscInt right)
 16: {
 17:   PetscInt  i,last;
 18:   PetscReal vl,tmp;

 21:   if (right <= 1) {
 22:     if (right == 1) {
 23:       if (v[0] > v[1]) SWAP(v[0],v[1],tmp);
 24:     }
 25:     return(0);
 26:   }
 27:   SWAP(v[0],v[right/2],tmp);
 28:   vl   = v[0];
 29:   last = 0;
 30:   for (i=1; i<=right; i++) {
 31:     if (v[i] < vl) {last++; SWAP(v[last],v[i],tmp);}
 32:   }
 33:   SWAP(v[0],v[last],tmp);
 34:   PetscSortReal_Private(v,last-1);
 35:   PetscSortReal_Private(v+last+1,right-(last+1));
 36:   return(0);
 37: }

 41: /*@
 42:    PetscSortReal - Sorts an array of doubles in place in increasing order.

 44:    Not Collective

 46:    Input Parameters:
 47: +  n  - number of values
 48: -  v  - array of doubles

 50:    Level: intermediate

 52:    Concepts: sorting^doubles

 54: .seealso: PetscSortInt(), PetscSortRealWithPermutation()
 55: @*/
 56: PetscErrorCode  PetscSortReal(PetscInt n,PetscReal v[])
 57: {
 58:   PetscInt  j,k;
 59:   PetscReal tmp,vk;

 62:   if (n<8) {
 63:     for (k=0; k<n; k++) {
 64:       vk = v[k];
 65:       for (j=k+1; j<n; j++) {
 66:         if (vk > v[j]) {
 67:           SWAP(v[k],v[j],tmp);
 68:           vk = v[k];
 69:         }
 70:       }
 71:     }
 72:   } else PetscSortReal_Private(v,n-1);
 73:   return(0);
 74: }


 79: /*@
 80:    PetscSortRemoveDupsReal - Sorts an array of doubles in place in increasing order removes all duplicate entries

 82:    Not Collective

 84:    Input Parameters:
 85: +  n  - number of values
 86: -  v  - array of doubles

 88:    Output Parameter:
 89: .  n - number of non-redundant values

 91:    Level: intermediate

 93:    Concepts: sorting^doubles

 95: .seealso: PetscSortReal(), PetscSortRemoveDupsInt()
 96: @*/
 97: PetscErrorCode  PetscSortRemoveDupsReal(PetscInt *n,PetscReal v[])
 98: {
100:   PetscInt       i,s = 0,N = *n, b = 0;

103:   PetscSortReal(N,v);
104:   for (i=0; i<N-1; i++) {
105:     if (v[b+s+1] != v[b]) {
106:       v[b+1] = v[b+s+1]; b++;
107:     } else s++;
108:   }
109:   *n = N - s;
110:   return(0);
111: }

115: /*@
116:    PetscSortSplit - Quick-sort split of an array of PetscScalars in place.

118:    Not Collective

120:    Input Parameters:
121: +  ncut  - splitig index
122: .  n     - number of values to sort
123: .  a     - array of values
124: -  idx   - index for array a

126:    Output Parameters:
127: +  a     - permuted array of values such that its elements satisfy:
128:            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
129:            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
130: -  idx   - permuted index of array a

132:    Level: intermediate

134:    Concepts: sorting^doubles

136: .seealso: PetscSortInt(), PetscSortRealWithPermutation()
137: @*/
138: PetscErrorCode  PetscSortSplit(PetscInt ncut,PetscInt n,PetscScalar a[],PetscInt idx[])
139: {
140:   PetscInt    i,mid,last,itmp,j,first;
141:   PetscScalar d,tmp;
142:   PetscReal   abskey;

145:   first = 0;
146:   last  = n-1;
147:   if (ncut < first || ncut > last) return(0);

149:   while (1) {
150:     mid    = first;
151:     d      = a[mid];
152:     abskey = PetscAbsScalar(d);
153:     i      = last;
154:     for (j = first + 1; j <= i; ++j) {
155:       d = a[j];
156:       if (PetscAbsScalar(d) >= abskey) {
157:         ++mid;
158:         /* interchange */
159:         tmp = a[mid];  itmp = idx[mid];
160:         a[mid] = a[j]; idx[mid] = idx[j];
161:         a[j] = tmp;    idx[j] = itmp;
162:       }
163:     }

165:     /* interchange */
166:     tmp = a[mid];      itmp = idx[mid];
167:     a[mid] = a[first]; idx[mid] = idx[first];
168:     a[first] = tmp;    idx[first] = itmp;

170:     /* test for while loop */
171:     if (mid == ncut) break;
172:     else if (mid > ncut) last = mid - 1;
173:     else first = mid + 1;
174:   }
175:   return(0);
176: }

180: /*@
181:    PetscSortSplitReal - Quick-sort split of an array of PetscReals in place.

183:    Not Collective

185:    Input Parameters:
186: +  ncut  - splitig index
187: .  n     - number of values to sort
188: .  a     - array of values in PetscReal
189: -  idx   - index for array a

191:    Output Parameters:
192: +  a     - permuted array of real values such that its elements satisfy:
193:            abs(a[i]) >= abs(a[ncut-1]) for i < ncut and
194:            abs(a[i]) <= abs(a[ncut-1]) for i >= ncut
195: -  idx   - permuted index of array a

197:    Level: intermediate

199:    Concepts: sorting^doubles

201: .seealso: PetscSortInt(), PetscSortRealWithPermutation()
202: @*/
203: PetscErrorCode  PetscSortSplitReal(PetscInt ncut,PetscInt n,PetscReal a[],PetscInt idx[])
204: {
205:   PetscInt  i,mid,last,itmp,j,first;
206:   PetscReal d,tmp;
207:   PetscReal abskey;

210:   first = 0;
211:   last  = n-1;
212:   if (ncut < first || ncut > last) return(0);

214:   while (1) {
215:     mid    = first;
216:     d      = a[mid];
217:     abskey = PetscAbsReal(d);
218:     i      = last;
219:     for (j = first + 1; j <= i; ++j) {
220:       d = a[j];
221:       if (PetscAbsReal(d) >= abskey) {
222:         ++mid;
223:         /* interchange */
224:         tmp = a[mid];  itmp = idx[mid];
225:         a[mid] = a[j]; idx[mid] = idx[j];
226:         a[j] = tmp;    idx[j] = itmp;
227:       }
228:     }

230:     /* interchange */
231:     tmp = a[mid];      itmp = idx[mid];
232:     a[mid] = a[first]; idx[mid] = idx[first];
233:     a[first] = tmp;    idx[first] = itmp;

235:     /* test for while loop */
236:     if (mid == ncut) break;
237:     else if (mid > ncut) last = mid - 1;
238:     else first = mid + 1;
239:   }
240:   return(0);
241: }