1e5c89e4eSSatish Balay #define PETSC_DLL 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 */ 8e5c89e4eSSatish Balay #include "petsc.h" /*I "petsc.h" I*/ 9e5c89e4eSSatish Balay #include "petscsys.h" /*I "petscsys.h" I*/ 10e5c89e4eSSatish Balay 11e5c89e4eSSatish Balay #define SWAP(a,b,t) {t=a;a=b;b=t;} 12e5c89e4eSSatish Balay 13e5c89e4eSSatish Balay #undef __FUNCT__ 14e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortReal_Private" 15e5c89e4eSSatish Balay /* A simple version of quicksort; taken from Kernighan and Ritchie, page 87 */ 16e5c89e4eSSatish Balay static PetscErrorCode PetscSortReal_Private(PetscReal *v,PetscInt right) 17e5c89e4eSSatish Balay { 18e5c89e4eSSatish Balay PetscInt i,last; 19e5c89e4eSSatish Balay PetscReal vl,tmp; 20e5c89e4eSSatish Balay 21e5c89e4eSSatish Balay PetscFunctionBegin; 22e5c89e4eSSatish Balay if (right <= 1) { 23e5c89e4eSSatish Balay if (right == 1) { 24e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP(v[0],v[1],tmp); 25e5c89e4eSSatish Balay } 26e5c89e4eSSatish Balay PetscFunctionReturn(0); 27e5c89e4eSSatish Balay } 28e5c89e4eSSatish Balay SWAP(v[0],v[right/2],tmp); 29e5c89e4eSSatish Balay vl = v[0]; 30e5c89e4eSSatish Balay last = 0; 31e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 32e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP(v[last],v[i],tmp);} 33e5c89e4eSSatish Balay } 34e5c89e4eSSatish Balay SWAP(v[0],v[last],tmp); 35e5c89e4eSSatish Balay PetscSortReal_Private(v,last-1); 36e5c89e4eSSatish Balay PetscSortReal_Private(v+last+1,right-(last+1)); 37e5c89e4eSSatish Balay PetscFunctionReturn(0); 38e5c89e4eSSatish Balay } 39e5c89e4eSSatish Balay 40e5c89e4eSSatish Balay #undef __FUNCT__ 41e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortReal" 42e5c89e4eSSatish Balay /*@ 43e5c89e4eSSatish Balay PetscSortReal - Sorts an array of doubles in place in increasing order. 44e5c89e4eSSatish Balay 45e5c89e4eSSatish Balay Not Collective 46e5c89e4eSSatish Balay 47e5c89e4eSSatish Balay Input Parameters: 48e5c89e4eSSatish Balay + n - number of values 49e5c89e4eSSatish Balay - v - array of doubles 50e5c89e4eSSatish Balay 51e5c89e4eSSatish Balay Level: intermediate 52e5c89e4eSSatish Balay 53e5c89e4eSSatish Balay Concepts: sorting^doubles 54e5c89e4eSSatish Balay 55e5c89e4eSSatish Balay .seealso: PetscSortInt(), PetscSortRealWithPermutation() 56e5c89e4eSSatish Balay @*/ 57e5c89e4eSSatish Balay PetscErrorCode PETSC_DLLEXPORT PetscSortReal(PetscInt n,PetscReal v[]) 58e5c89e4eSSatish Balay { 59e5c89e4eSSatish Balay PetscInt j,k; 60e5c89e4eSSatish Balay PetscReal tmp,vk; 61e5c89e4eSSatish Balay 62e5c89e4eSSatish Balay PetscFunctionBegin; 63e5c89e4eSSatish Balay if (n<8) { 64e5c89e4eSSatish Balay for (k=0; k<n; k++) { 65e5c89e4eSSatish Balay vk = v[k]; 66e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 67e5c89e4eSSatish Balay if (vk > v[j]) { 68e5c89e4eSSatish Balay SWAP(v[k],v[j],tmp); 69e5c89e4eSSatish Balay vk = v[k]; 70e5c89e4eSSatish Balay } 71e5c89e4eSSatish Balay } 72e5c89e4eSSatish Balay } 73e5c89e4eSSatish Balay } else { 74e5c89e4eSSatish Balay PetscSortReal_Private(v,n-1); 75e5c89e4eSSatish Balay } 76e5c89e4eSSatish Balay PetscFunctionReturn(0); 77e5c89e4eSSatish Balay } 78e5c89e4eSSatish Balay 79d97c5584SHong Zhang #undef __FUNCT__ 80d97c5584SHong Zhang #define __FUNCT__ "PetscSortSplit" 81d97c5584SHong Zhang /*@ 82*021822bcSHong Zhang PetscSortSplit - Quick-sort split of an array of PetscScalars in place. 83d97c5584SHong Zhang 84d97c5584SHong Zhang Not Collective 85d97c5584SHong Zhang 86d97c5584SHong Zhang Input Parameters: 87d97c5584SHong Zhang + ncut - splitig index 88d97c5584SHong Zhang . n - number of values to sort 89d97c5584SHong Zhang . a - array of values 90d97c5584SHong Zhang - idx - index for array a 91d97c5584SHong Zhang 92d97c5584SHong Zhang Output Parameters: 93d97c5584SHong Zhang + a - permuted array of values such that its elements satisfy: 94d97c5584SHong Zhang abs(a[i]) >= abs(a[ncut-1]) for i < ncut and 95d97c5584SHong Zhang abs(a[i]) <= abs(a[ncut-1]) for i >= ncut 96d97c5584SHong Zhang - idx - permuted index of array a 97d97c5584SHong Zhang 98d97c5584SHong Zhang Level: intermediate 99d97c5584SHong Zhang 100d97c5584SHong Zhang Concepts: sorting^doubles 101d97c5584SHong Zhang 102d97c5584SHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation() 103d97c5584SHong Zhang @*/ 104d97c5584SHong Zhang PetscErrorCode PETSC_DLLEXPORT PetscSortSplit(PetscInt ncut,PetscInt n,PetscScalar a[],PetscInt idx[]) 105d97c5584SHong Zhang { 106d97c5584SHong Zhang PetscInt i,mid,last,itmp,j,first; 107d97c5584SHong Zhang PetscScalar d,tmp; 108d97c5584SHong Zhang PetscReal abskey; 109d97c5584SHong Zhang 110d97c5584SHong Zhang PetscFunctionBegin; 111d97c5584SHong Zhang first = 0; 112d97c5584SHong Zhang last = n-1; 113d97c5584SHong Zhang if (ncut < first || ncut > last) PetscFunctionReturn(0); 114d97c5584SHong Zhang 115d97c5584SHong Zhang while (1){ 116d97c5584SHong Zhang mid = first; 117d97c5584SHong Zhang abskey = (d = a[mid],PetscAbsScalar(d)); 118d97c5584SHong Zhang i = last; 119d97c5584SHong Zhang for (j = first + 1; j <= i; ++j) { 120d97c5584SHong Zhang if ((d = a[j],PetscAbsScalar(d)) >= abskey) { 121d97c5584SHong Zhang ++mid; 122d97c5584SHong Zhang /* interchange */ 123d97c5584SHong Zhang tmp = a[mid]; itmp = idx[mid]; 124d97c5584SHong Zhang a[mid] = a[j]; idx[mid] = idx[j]; 125d97c5584SHong Zhang a[j] = tmp; idx[j] = itmp; 126d97c5584SHong Zhang } 127d97c5584SHong Zhang } 128d97c5584SHong Zhang 129d97c5584SHong Zhang /* interchange */ 130d97c5584SHong Zhang tmp = a[mid]; itmp = idx[mid]; 131d97c5584SHong Zhang a[mid] = a[first]; idx[mid] = idx[first]; 132d97c5584SHong Zhang a[first] = tmp; idx[first] = itmp; 133d97c5584SHong Zhang 134d97c5584SHong Zhang /* test for while loop */ 135d97c5584SHong Zhang if (mid == ncut) { 136d97c5584SHong Zhang break; 137d97c5584SHong Zhang } else if (mid > ncut){ 138d97c5584SHong Zhang last = mid - 1; 139d97c5584SHong Zhang } else { 140d97c5584SHong Zhang first = mid + 1; 141d97c5584SHong Zhang } 142d97c5584SHong Zhang } 143d97c5584SHong Zhang PetscFunctionReturn(0); 144d97c5584SHong Zhang } 145*021822bcSHong Zhang 146*021822bcSHong Zhang #undef __FUNCT__ 147*021822bcSHong Zhang #define __FUNCT__ "PetscSortSplitReal" 148*021822bcSHong Zhang /*@ 149*021822bcSHong Zhang PetscSortSplitReal - Quick-sort split of an array of PetscReals in place. 150*021822bcSHong Zhang 151*021822bcSHong Zhang Not Collective 152*021822bcSHong Zhang 153*021822bcSHong Zhang Input Parameters: 154*021822bcSHong Zhang + ncut - splitig index 155*021822bcSHong Zhang . n - number of values to sort 156*021822bcSHong Zhang . a - array of values in PetscReal 157*021822bcSHong Zhang - idx - index for array a 158*021822bcSHong Zhang 159*021822bcSHong Zhang Output Parameters: 160*021822bcSHong Zhang + a - permuted array of real values such that its elements satisfy: 161*021822bcSHong Zhang abs(a[i]) >= abs(a[ncut-1]) for i < ncut and 162*021822bcSHong Zhang abs(a[i]) <= abs(a[ncut-1]) for i >= ncut 163*021822bcSHong Zhang - idx - permuted index of array a 164*021822bcSHong Zhang 165*021822bcSHong Zhang Level: intermediate 166*021822bcSHong Zhang 167*021822bcSHong Zhang Concepts: sorting^doubles 168*021822bcSHong Zhang 169*021822bcSHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation() 170*021822bcSHong Zhang @*/ 171*021822bcSHong Zhang PetscErrorCode PETSC_DLLEXPORT PetscSortSplitReal(PetscInt ncut,PetscInt n,PetscReal a[],PetscInt idx[]) 172*021822bcSHong Zhang { 173*021822bcSHong Zhang PetscInt i,mid,last,itmp,j,first; 174*021822bcSHong Zhang PetscReal d,tmp; 175*021822bcSHong Zhang PetscReal abskey; 176*021822bcSHong Zhang 177*021822bcSHong Zhang PetscFunctionBegin; 178*021822bcSHong Zhang first = 0; 179*021822bcSHong Zhang last = n-1; 180*021822bcSHong Zhang if (ncut < first || ncut > last) PetscFunctionReturn(0); 181*021822bcSHong Zhang 182*021822bcSHong Zhang while (1){ 183*021822bcSHong Zhang mid = first; 184*021822bcSHong Zhang abskey = (d = a[mid],PetscAbsScalar(d)); 185*021822bcSHong Zhang i = last; 186*021822bcSHong Zhang for (j = first + 1; j <= i; ++j) { 187*021822bcSHong Zhang if ((d = a[j],PetscAbsScalar(d)) >= abskey) { 188*021822bcSHong Zhang ++mid; 189*021822bcSHong Zhang /* interchange */ 190*021822bcSHong Zhang tmp = a[mid]; itmp = idx[mid]; 191*021822bcSHong Zhang a[mid] = a[j]; idx[mid] = idx[j]; 192*021822bcSHong Zhang a[j] = tmp; idx[j] = itmp; 193*021822bcSHong Zhang } 194*021822bcSHong Zhang } 195*021822bcSHong Zhang 196*021822bcSHong Zhang /* interchange */ 197*021822bcSHong Zhang tmp = a[mid]; itmp = idx[mid]; 198*021822bcSHong Zhang a[mid] = a[first]; idx[mid] = idx[first]; 199*021822bcSHong Zhang a[first] = tmp; idx[first] = itmp; 200*021822bcSHong Zhang 201*021822bcSHong Zhang /* test for while loop */ 202*021822bcSHong Zhang if (mid == ncut) { 203*021822bcSHong Zhang break; 204*021822bcSHong Zhang } else if (mid > ncut){ 205*021822bcSHong Zhang last = mid - 1; 206*021822bcSHong Zhang } else { 207*021822bcSHong Zhang first = mid + 1; 208*021822bcSHong Zhang } 209*021822bcSHong Zhang } 210*021822bcSHong Zhang PetscFunctionReturn(0); 211*021822bcSHong Zhang } 212*021822bcSHong Zhang 213