Actual source code: sortd.c

  1: /*$Id: sortd.c,v 1.29 2001/08/07 21:29:06 bsmith Exp $*/
  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:    The word "register"  in this code is used to identify data that is not
  8:    aliased.  For some compilers, this can cause the compiler to fail to
  9:    place inner-loop variables into registers.
 10:  */
 11:  #include petsc.h
 12:  #include petscsys.h

 14: #define SWAP(a,b,t) {t=a;a=b;b=t;}
 15: 
 16: #undef __FUNCT__  
 18: /* A simple version of quicksort; taken from Kernighan and Ritchie, page 87 */
 19: static int PetsciDqsort(PetscReal *v,int right)
 20: {
 21:   int       i,last;
 22:   PetscReal vl,tmp;
 23: 
 25:   if (right <= 1) {
 26:     if (right == 1) {
 27:       if (v[0] > v[1]) SWAP(v[0],v[1],tmp);
 28:     }
 29:     return(0);
 30:   }
 31:   SWAP(v[0],v[right/2],tmp);
 32:   vl   = v[0];
 33:   last = 0;
 34:   for (i=1; i<=right; i++) {
 35:     if (v[i] < vl) {last++; SWAP(v[last],v[i],tmp);}
 36:   }
 37:   SWAP(v[0],v[last],tmp);
 38:   PetsciDqsort(v,last-1);
 39:   PetsciDqsort(v+last+1,right-(last+1));
 40:   return(0);
 41: }

 43: #undef __FUNCT__  
 45: /*@
 46:    PetscSortReal - Sorts an array of doubles in place in increasing order.

 48:    Not Collective

 50:    Input Parameters:
 51: +  n  - number of values
 52: -  v  - array of doubles

 54:    Level: intermediate

 56:    Concepts: sorting^doubles

 58: .seealso: PetscSortInt(), PetscSortRealWithPermutation()
 59: @*/
 60: int PetscSortReal(int n,PetscReal v[])
 61: {
 62:   int       j,k;
 63:   PetscReal tmp,vk;

 66:   if (n<8) {
 67:     for (k=0; k<n; k++) {
 68:       vk = v[k];
 69:       for (j=k+1; j<n; j++) {
 70:         if (vk > v[j]) {
 71:           SWAP(v[k],v[j],tmp);
 72:           vk = v[k];
 73:         }
 74:       }
 75:     }
 76:   } else {
 77:     PetsciDqsort(v,n-1);
 78:   }
 79:   return(0);
 80: }