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*/ 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 } 72a297a907SKarl Rupp } else PetscSortReal_Private(v,n-1); 73e5c89e4eSSatish Balay PetscFunctionReturn(0); 74e5c89e4eSSatish Balay } 75e5c89e4eSSatish Balay 76745b41b2SMatthew G. Knepley 77745b41b2SMatthew G. Knepley #undef __FUNCT__ 78745b41b2SMatthew G. Knepley #define __FUNCT__ "PetscSortRemoveDupsReal" 79745b41b2SMatthew G. Knepley /*@ 80745b41b2SMatthew G. Knepley PetscSortRemoveDupsReal - Sorts an array of doubles in place in increasing order removes all duplicate entries 81745b41b2SMatthew G. Knepley 82745b41b2SMatthew G. Knepley Not Collective 83745b41b2SMatthew G. Knepley 84745b41b2SMatthew G. Knepley Input Parameters: 85745b41b2SMatthew G. Knepley + n - number of values 86745b41b2SMatthew G. Knepley - v - array of doubles 87745b41b2SMatthew G. Knepley 88745b41b2SMatthew G. Knepley Output Parameter: 89745b41b2SMatthew G. Knepley . n - number of non-redundant values 90745b41b2SMatthew G. Knepley 91745b41b2SMatthew G. Knepley Level: intermediate 92745b41b2SMatthew G. Knepley 93745b41b2SMatthew G. Knepley Concepts: sorting^doubles 94745b41b2SMatthew G. Knepley 95745b41b2SMatthew G. Knepley .seealso: PetscSortReal(), PetscSortRemoveDupsInt() 96745b41b2SMatthew G. Knepley @*/ 97745b41b2SMatthew G. Knepley PetscErrorCode PetscSortRemoveDupsReal(PetscInt *n,PetscReal v[]) 98745b41b2SMatthew G. Knepley { 99745b41b2SMatthew G. Knepley PetscErrorCode ierr; 100745b41b2SMatthew G. Knepley PetscInt i,s = 0,N = *n, b = 0; 101745b41b2SMatthew G. Knepley 102745b41b2SMatthew G. Knepley PetscFunctionBegin; 103745b41b2SMatthew G. Knepley ierr = PetscSortReal(N,v);CHKERRQ(ierr); 104745b41b2SMatthew G. Knepley for (i=0; i<N-1; i++) { 105745b41b2SMatthew G. Knepley if (v[b+s+1] != v[b]) { 106745b41b2SMatthew G. Knepley v[b+1] = v[b+s+1]; b++; 107745b41b2SMatthew G. Knepley } else s++; 108745b41b2SMatthew G. Knepley } 109745b41b2SMatthew G. Knepley *n = N - s; 110745b41b2SMatthew G. Knepley PetscFunctionReturn(0); 111745b41b2SMatthew G. Knepley } 112745b41b2SMatthew G. Knepley 113d97c5584SHong Zhang #undef __FUNCT__ 114d97c5584SHong Zhang #define __FUNCT__ "PetscSortSplit" 115d97c5584SHong Zhang /*@ 116021822bcSHong Zhang PetscSortSplit - Quick-sort split of an array of PetscScalars in place. 117d97c5584SHong Zhang 118d97c5584SHong Zhang Not Collective 119d97c5584SHong Zhang 120d97c5584SHong Zhang Input Parameters: 121d97c5584SHong Zhang + ncut - splitig index 122d97c5584SHong Zhang . n - number of values to sort 123d97c5584SHong Zhang . a - array of values 124d97c5584SHong Zhang - idx - index for array a 125d97c5584SHong Zhang 126d97c5584SHong Zhang Output Parameters: 127d97c5584SHong Zhang + a - permuted array of values such that its elements satisfy: 128d97c5584SHong Zhang abs(a[i]) >= abs(a[ncut-1]) for i < ncut and 129d97c5584SHong Zhang abs(a[i]) <= abs(a[ncut-1]) for i >= ncut 130d97c5584SHong Zhang - idx - permuted index of array a 131d97c5584SHong Zhang 132d97c5584SHong Zhang Level: intermediate 133d97c5584SHong Zhang 134d97c5584SHong Zhang Concepts: sorting^doubles 135d97c5584SHong Zhang 136d97c5584SHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation() 137d97c5584SHong Zhang @*/ 1387087cfbeSBarry Smith PetscErrorCode PetscSortSplit(PetscInt ncut,PetscInt n,PetscScalar a[],PetscInt idx[]) 139d97c5584SHong Zhang { 140d97c5584SHong Zhang PetscInt i,mid,last,itmp,j,first; 141d97c5584SHong Zhang PetscScalar d,tmp; 142d97c5584SHong Zhang PetscReal abskey; 143d97c5584SHong Zhang 144d97c5584SHong Zhang PetscFunctionBegin; 145d97c5584SHong Zhang first = 0; 146d97c5584SHong Zhang last = n-1; 147d97c5584SHong Zhang if (ncut < first || ncut > last) PetscFunctionReturn(0); 148d97c5584SHong Zhang 149d97c5584SHong Zhang while (1) { 150d97c5584SHong Zhang mid = first; 151*2a274a78SSatish Balay d = a[mid]; 152*2a274a78SSatish Balay abskey = PetscAbsScalar(d); 153d97c5584SHong Zhang i = last; 154d97c5584SHong Zhang for (j = first + 1; j <= i; ++j) { 155*2a274a78SSatish Balay d = a[j]; 156*2a274a78SSatish Balay if (PetscAbsScalar(d) >= abskey) { 157d97c5584SHong Zhang ++mid; 158d97c5584SHong Zhang /* interchange */ 159d97c5584SHong Zhang tmp = a[mid]; itmp = idx[mid]; 160d97c5584SHong Zhang a[mid] = a[j]; idx[mid] = idx[j]; 161d97c5584SHong Zhang a[j] = tmp; idx[j] = itmp; 162d97c5584SHong Zhang } 163d97c5584SHong Zhang } 164d97c5584SHong Zhang 165d97c5584SHong Zhang /* interchange */ 166d97c5584SHong Zhang tmp = a[mid]; itmp = idx[mid]; 167d97c5584SHong Zhang a[mid] = a[first]; idx[mid] = idx[first]; 168d97c5584SHong Zhang a[first] = tmp; idx[first] = itmp; 169d97c5584SHong Zhang 170d97c5584SHong Zhang /* test for while loop */ 171a297a907SKarl Rupp if (mid == ncut) break; 172a297a907SKarl Rupp else if (mid > ncut) last = mid - 1; 173a297a907SKarl Rupp else first = mid + 1; 174d97c5584SHong Zhang } 175d97c5584SHong Zhang PetscFunctionReturn(0); 176d97c5584SHong Zhang } 177021822bcSHong Zhang 178021822bcSHong Zhang #undef __FUNCT__ 179021822bcSHong Zhang #define __FUNCT__ "PetscSortSplitReal" 180021822bcSHong Zhang /*@ 181021822bcSHong Zhang PetscSortSplitReal - Quick-sort split of an array of PetscReals in place. 182021822bcSHong Zhang 183021822bcSHong Zhang Not Collective 184021822bcSHong Zhang 185021822bcSHong Zhang Input Parameters: 186021822bcSHong Zhang + ncut - splitig index 187021822bcSHong Zhang . n - number of values to sort 188021822bcSHong Zhang . a - array of values in PetscReal 189021822bcSHong Zhang - idx - index for array a 190021822bcSHong Zhang 191021822bcSHong Zhang Output Parameters: 192021822bcSHong Zhang + a - permuted array of real values such that its elements satisfy: 193021822bcSHong Zhang abs(a[i]) >= abs(a[ncut-1]) for i < ncut and 194021822bcSHong Zhang abs(a[i]) <= abs(a[ncut-1]) for i >= ncut 195021822bcSHong Zhang - idx - permuted index of array a 196021822bcSHong Zhang 197021822bcSHong Zhang Level: intermediate 198021822bcSHong Zhang 199021822bcSHong Zhang Concepts: sorting^doubles 200021822bcSHong Zhang 201021822bcSHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation() 202021822bcSHong Zhang @*/ 2037087cfbeSBarry Smith PetscErrorCode PetscSortSplitReal(PetscInt ncut,PetscInt n,PetscReal a[],PetscInt idx[]) 204021822bcSHong Zhang { 205021822bcSHong Zhang PetscInt i,mid,last,itmp,j,first; 206021822bcSHong Zhang PetscReal d,tmp; 207021822bcSHong Zhang PetscReal abskey; 208021822bcSHong Zhang 209021822bcSHong Zhang PetscFunctionBegin; 210021822bcSHong Zhang first = 0; 211021822bcSHong Zhang last = n-1; 212021822bcSHong Zhang if (ncut < first || ncut > last) PetscFunctionReturn(0); 213021822bcSHong Zhang 214021822bcSHong Zhang while (1) { 215021822bcSHong Zhang mid = first; 216*2a274a78SSatish Balay d = a[mid]; 217*2a274a78SSatish Balay abskey = PetscAbsReal(d); 218021822bcSHong Zhang i = last; 219021822bcSHong Zhang for (j = first + 1; j <= i; ++j) { 220*2a274a78SSatish Balay d = a[j]; 221*2a274a78SSatish Balay if (PetscAbsReal(d) >= abskey) { 222021822bcSHong Zhang ++mid; 223021822bcSHong Zhang /* interchange */ 224021822bcSHong Zhang tmp = a[mid]; itmp = idx[mid]; 225021822bcSHong Zhang a[mid] = a[j]; idx[mid] = idx[j]; 226021822bcSHong Zhang a[j] = tmp; idx[j] = itmp; 227021822bcSHong Zhang } 228021822bcSHong Zhang } 229021822bcSHong Zhang 230021822bcSHong Zhang /* interchange */ 231021822bcSHong Zhang tmp = a[mid]; itmp = idx[mid]; 232021822bcSHong Zhang a[mid] = a[first]; idx[mid] = idx[first]; 233021822bcSHong Zhang a[first] = tmp; idx[first] = itmp; 234021822bcSHong Zhang 235021822bcSHong Zhang /* test for while loop */ 236a297a907SKarl Rupp if (mid == ncut) break; 237a297a907SKarl Rupp else if (mid > ncut) last = mid - 1; 238a297a907SKarl Rupp else first = mid + 1; 239021822bcSHong Zhang } 240021822bcSHong Zhang PetscFunctionReturn(0); 241021822bcSHong Zhang } 242021822bcSHong Zhang 243