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*/ 939f41f7fSStefano Zampini #include <petsc/private/petscimpl.h> 10e5c89e4eSSatish Balay 11e5c89e4eSSatish Balay #define SWAP(a,b,t) {t=a;a=b;b=t;} 12e5c89e4eSSatish Balay 13*6a7c706bSVaclav Hapla /*@ 14*6a7c706bSVaclav Hapla PetscSortedReal - Determines whether the array is sorted. 15*6a7c706bSVaclav Hapla 16*6a7c706bSVaclav Hapla Not Collective 17*6a7c706bSVaclav Hapla 18*6a7c706bSVaclav Hapla Input Parameters: 19*6a7c706bSVaclav Hapla + n - number of values 20*6a7c706bSVaclav Hapla - X - array of integers 21*6a7c706bSVaclav Hapla 22*6a7c706bSVaclav Hapla Output Parameters: 23*6a7c706bSVaclav Hapla . sorted - flag whether the array is sorted 24*6a7c706bSVaclav Hapla 25*6a7c706bSVaclav Hapla Level: intermediate 26*6a7c706bSVaclav Hapla 27*6a7c706bSVaclav Hapla .seealso: PetscSortReal(), PetscSortedInt(), PetscSortedMPIInt() 28*6a7c706bSVaclav Hapla @*/ 29*6a7c706bSVaclav Hapla PetscErrorCode PetscSortedReal(PetscInt n,const PetscReal X[],PetscBool *sorted) 30*6a7c706bSVaclav Hapla { 31*6a7c706bSVaclav Hapla PetscFunctionBegin; 32*6a7c706bSVaclav Hapla PetscSorted(n,X,*sorted); 33*6a7c706bSVaclav Hapla PetscFunctionReturn(0); 34*6a7c706bSVaclav Hapla } 35*6a7c706bSVaclav Hapla 36e5c89e4eSSatish Balay /* A simple version of quicksort; taken from Kernighan and Ritchie, page 87 */ 37e5c89e4eSSatish Balay static PetscErrorCode PetscSortReal_Private(PetscReal *v,PetscInt right) 38e5c89e4eSSatish Balay { 39e5c89e4eSSatish Balay PetscInt i,last; 40e5c89e4eSSatish Balay PetscReal vl,tmp; 41e5c89e4eSSatish Balay 42e5c89e4eSSatish Balay PetscFunctionBegin; 43e5c89e4eSSatish Balay if (right <= 1) { 44e5c89e4eSSatish Balay if (right == 1) { 45e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP(v[0],v[1],tmp); 46e5c89e4eSSatish Balay } 47e5c89e4eSSatish Balay PetscFunctionReturn(0); 48e5c89e4eSSatish Balay } 49e5c89e4eSSatish Balay SWAP(v[0],v[right/2],tmp); 50e5c89e4eSSatish Balay vl = v[0]; 51e5c89e4eSSatish Balay last = 0; 52e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 53e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP(v[last],v[i],tmp);} 54e5c89e4eSSatish Balay } 55e5c89e4eSSatish Balay SWAP(v[0],v[last],tmp); 56e5c89e4eSSatish Balay PetscSortReal_Private(v,last-1); 57e5c89e4eSSatish Balay PetscSortReal_Private(v+last+1,right-(last+1)); 58e5c89e4eSSatish Balay PetscFunctionReturn(0); 59e5c89e4eSSatish Balay } 60e5c89e4eSSatish Balay 61e5c89e4eSSatish Balay /*@ 62e5c89e4eSSatish Balay PetscSortReal - Sorts an array of doubles in place in increasing order. 63e5c89e4eSSatish Balay 64e5c89e4eSSatish Balay Not Collective 65e5c89e4eSSatish Balay 66e5c89e4eSSatish Balay Input Parameters: 67e5c89e4eSSatish Balay + n - number of values 68e5c89e4eSSatish Balay - v - array of doubles 69e5c89e4eSSatish Balay 70e5c89e4eSSatish Balay Level: intermediate 71e5c89e4eSSatish Balay 7239f41f7fSStefano Zampini .seealso: PetscSortInt(), PetscSortRealWithPermutation(), PetscSortRealWithArrayInt() 73e5c89e4eSSatish Balay @*/ 747087cfbeSBarry Smith PetscErrorCode PetscSortReal(PetscInt n,PetscReal v[]) 75e5c89e4eSSatish Balay { 76e5c89e4eSSatish Balay PetscInt j,k; 77e5c89e4eSSatish Balay PetscReal tmp,vk; 78e5c89e4eSSatish Balay 79e5c89e4eSSatish Balay PetscFunctionBegin; 8039f41f7fSStefano Zampini PetscValidPointer(v,2); 81e5c89e4eSSatish Balay if (n<8) { 82e5c89e4eSSatish Balay for (k=0; k<n; k++) { 83e5c89e4eSSatish Balay vk = v[k]; 84e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 85e5c89e4eSSatish Balay if (vk > v[j]) { 86e5c89e4eSSatish Balay SWAP(v[k],v[j],tmp); 87e5c89e4eSSatish Balay vk = v[k]; 88e5c89e4eSSatish Balay } 89e5c89e4eSSatish Balay } 90e5c89e4eSSatish Balay } 91a297a907SKarl Rupp } else PetscSortReal_Private(v,n-1); 92e5c89e4eSSatish Balay PetscFunctionReturn(0); 93e5c89e4eSSatish Balay } 94e5c89e4eSSatish Balay 9539f41f7fSStefano Zampini #define SWAP2ri(a,b,c,d,rt,it) {rt=a;a=b;b=rt;it=c;c=d;d=it;} 9639f41f7fSStefano Zampini 9739f41f7fSStefano Zampini /* modified from PetscSortIntWithArray_Private */ 9839f41f7fSStefano Zampini static PetscErrorCode PetscSortRealWithArrayInt_Private(PetscReal *v,PetscInt *V,PetscInt right) 9939f41f7fSStefano Zampini { 10039f41f7fSStefano Zampini PetscErrorCode ierr; 10139f41f7fSStefano Zampini PetscInt i,last,itmp; 10239f41f7fSStefano Zampini PetscReal rvl,rtmp; 10339f41f7fSStefano Zampini 10439f41f7fSStefano Zampini PetscFunctionBegin; 10539f41f7fSStefano Zampini if (right <= 1) { 10639f41f7fSStefano Zampini if (right == 1) { 10739f41f7fSStefano Zampini if (v[0] > v[1]) SWAP2ri(v[0],v[1],V[0],V[1],rtmp,itmp); 10839f41f7fSStefano Zampini } 10939f41f7fSStefano Zampini PetscFunctionReturn(0); 11039f41f7fSStefano Zampini } 11139f41f7fSStefano Zampini SWAP2ri(v[0],v[right/2],V[0],V[right/2],rtmp,itmp); 11239f41f7fSStefano Zampini rvl = v[0]; 11339f41f7fSStefano Zampini last = 0; 11439f41f7fSStefano Zampini for (i=1; i<=right; i++) { 11539f41f7fSStefano Zampini if (v[i] < rvl) {last++; SWAP2ri(v[last],v[i],V[last],V[i],rtmp,itmp);} 11639f41f7fSStefano Zampini } 11739f41f7fSStefano Zampini SWAP2ri(v[0],v[last],V[0],V[last],rtmp,itmp); 11839f41f7fSStefano Zampini ierr = PetscSortRealWithArrayInt_Private(v,V,last-1);CHKERRQ(ierr); 11939f41f7fSStefano Zampini ierr = PetscSortRealWithArrayInt_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 12039f41f7fSStefano Zampini PetscFunctionReturn(0); 12139f41f7fSStefano Zampini } 12239f41f7fSStefano Zampini /*@ 12339f41f7fSStefano Zampini PetscSortRealWithArrayInt - Sorts an array of PetscReal in place in increasing order; 12439f41f7fSStefano Zampini changes a second PetscInt array to match the sorted first array. 12539f41f7fSStefano Zampini 12639f41f7fSStefano Zampini Not Collective 12739f41f7fSStefano Zampini 12839f41f7fSStefano Zampini Input Parameters: 12939f41f7fSStefano Zampini + n - number of values 13039f41f7fSStefano Zampini . i - array of integers 13139f41f7fSStefano Zampini - I - second array of integers 13239f41f7fSStefano Zampini 13339f41f7fSStefano Zampini Level: intermediate 13439f41f7fSStefano Zampini 13539f41f7fSStefano Zampini .seealso: PetscSortReal() 13639f41f7fSStefano Zampini @*/ 13739f41f7fSStefano Zampini PetscErrorCode PetscSortRealWithArrayInt(PetscInt n,PetscReal r[],PetscInt Ii[]) 13839f41f7fSStefano Zampini { 13939f41f7fSStefano Zampini PetscErrorCode ierr; 14039f41f7fSStefano Zampini PetscInt j,k,itmp; 14139f41f7fSStefano Zampini PetscReal rk,rtmp; 14239f41f7fSStefano Zampini 14339f41f7fSStefano Zampini PetscFunctionBegin; 14439f41f7fSStefano Zampini PetscValidPointer(r,2); 14539f41f7fSStefano Zampini PetscValidPointer(Ii,3); 14639f41f7fSStefano Zampini if (n<8) { 14739f41f7fSStefano Zampini for (k=0; k<n; k++) { 14839f41f7fSStefano Zampini rk = r[k]; 14939f41f7fSStefano Zampini for (j=k+1; j<n; j++) { 15039f41f7fSStefano Zampini if (rk > r[j]) { 15139f41f7fSStefano Zampini SWAP2ri(r[k],r[j],Ii[k],Ii[j],rtmp,itmp); 15239f41f7fSStefano Zampini rk = r[k]; 15339f41f7fSStefano Zampini } 15439f41f7fSStefano Zampini } 15539f41f7fSStefano Zampini } 15639f41f7fSStefano Zampini } else { 15739f41f7fSStefano Zampini ierr = PetscSortRealWithArrayInt_Private(r,Ii,n-1);CHKERRQ(ierr); 15839f41f7fSStefano Zampini } 15939f41f7fSStefano Zampini PetscFunctionReturn(0); 16039f41f7fSStefano Zampini } 16139f41f7fSStefano Zampini 16239f41f7fSStefano Zampini /*@ 16339f41f7fSStefano Zampini PetscFindReal - Finds a PetscReal in a sorted array of PetscReals 16439f41f7fSStefano Zampini 16539f41f7fSStefano Zampini Not Collective 16639f41f7fSStefano Zampini 16739f41f7fSStefano Zampini Input Parameters: 16839f41f7fSStefano Zampini + key - the value to locate 16939f41f7fSStefano Zampini . n - number of values in the array 17039f41f7fSStefano Zampini . ii - array of values 17139f41f7fSStefano Zampini - eps - tolerance used to compare 17239f41f7fSStefano Zampini 17339f41f7fSStefano Zampini Output Parameter: 17439f41f7fSStefano Zampini . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 17539f41f7fSStefano Zampini 17639f41f7fSStefano Zampini Level: intermediate 17739f41f7fSStefano Zampini 17839f41f7fSStefano Zampini .seealso: PetscSortReal(), PetscSortRealWithArrayInt() 17939f41f7fSStefano Zampini @*/ 18039f41f7fSStefano Zampini PetscErrorCode PetscFindReal(PetscReal key, PetscInt n, const PetscReal t[], PetscReal eps, PetscInt *loc) 18139f41f7fSStefano Zampini { 18239f41f7fSStefano Zampini PetscInt lo = 0,hi = n; 18339f41f7fSStefano Zampini 18439f41f7fSStefano Zampini PetscFunctionBegin; 18539f41f7fSStefano Zampini PetscValidPointer(loc,4); 18639f41f7fSStefano Zampini if (!n) {*loc = -1; PetscFunctionReturn(0);} 18739f41f7fSStefano Zampini PetscValidPointer(t,3); 188*6a7c706bSVaclav Hapla PetscCheckSorted(n,t); 18939f41f7fSStefano Zampini while (hi - lo > 1) { 19039f41f7fSStefano Zampini PetscInt mid = lo + (hi - lo)/2; 19139f41f7fSStefano Zampini if (key < t[mid]) hi = mid; 19239f41f7fSStefano Zampini else lo = mid; 19339f41f7fSStefano Zampini } 19439f41f7fSStefano Zampini *loc = (PetscAbsReal(key - t[lo]) < eps) ? lo : -(lo + (key > t[lo]) + 1); 19539f41f7fSStefano Zampini PetscFunctionReturn(0); 19639f41f7fSStefano Zampini } 197745b41b2SMatthew G. Knepley 198745b41b2SMatthew G. Knepley /*@ 199745b41b2SMatthew G. Knepley PetscSortRemoveDupsReal - Sorts an array of doubles in place in increasing order removes all duplicate entries 200745b41b2SMatthew G. Knepley 201745b41b2SMatthew G. Knepley Not Collective 202745b41b2SMatthew G. Knepley 203745b41b2SMatthew G. Knepley Input Parameters: 204745b41b2SMatthew G. Knepley + n - number of values 205745b41b2SMatthew G. Knepley - v - array of doubles 206745b41b2SMatthew G. Knepley 207745b41b2SMatthew G. Knepley Output Parameter: 208745b41b2SMatthew G. Knepley . n - number of non-redundant values 209745b41b2SMatthew G. Knepley 210745b41b2SMatthew G. Knepley Level: intermediate 211745b41b2SMatthew G. Knepley 212745b41b2SMatthew G. Knepley .seealso: PetscSortReal(), PetscSortRemoveDupsInt() 213745b41b2SMatthew G. Knepley @*/ 214745b41b2SMatthew G. Knepley PetscErrorCode PetscSortRemoveDupsReal(PetscInt *n,PetscReal v[]) 215745b41b2SMatthew G. Knepley { 216745b41b2SMatthew G. Knepley PetscErrorCode ierr; 217745b41b2SMatthew G. Knepley PetscInt i,s = 0,N = *n, b = 0; 218745b41b2SMatthew G. Knepley 219745b41b2SMatthew G. Knepley PetscFunctionBegin; 220745b41b2SMatthew G. Knepley ierr = PetscSortReal(N,v);CHKERRQ(ierr); 221745b41b2SMatthew G. Knepley for (i=0; i<N-1; i++) { 222745b41b2SMatthew G. Knepley if (v[b+s+1] != v[b]) { 223745b41b2SMatthew G. Knepley v[b+1] = v[b+s+1]; b++; 224745b41b2SMatthew G. Knepley } else s++; 225745b41b2SMatthew G. Knepley } 226745b41b2SMatthew G. Knepley *n = N - s; 227745b41b2SMatthew G. Knepley PetscFunctionReturn(0); 228745b41b2SMatthew G. Knepley } 229745b41b2SMatthew G. Knepley 230d97c5584SHong Zhang /*@ 231021822bcSHong Zhang PetscSortSplit - Quick-sort split of an array of PetscScalars in place. 232d97c5584SHong Zhang 233d97c5584SHong Zhang Not Collective 234d97c5584SHong Zhang 235d97c5584SHong Zhang Input Parameters: 236d97c5584SHong Zhang + ncut - splitig index 237d97c5584SHong Zhang . n - number of values to sort 238d97c5584SHong Zhang . a - array of values 239d97c5584SHong Zhang - idx - index for array a 240d97c5584SHong Zhang 241d97c5584SHong Zhang Output Parameters: 242d97c5584SHong Zhang + a - permuted array of values such that its elements satisfy: 243d97c5584SHong Zhang abs(a[i]) >= abs(a[ncut-1]) for i < ncut and 244d97c5584SHong Zhang abs(a[i]) <= abs(a[ncut-1]) for i >= ncut 245d97c5584SHong Zhang - idx - permuted index of array a 246d97c5584SHong Zhang 247d97c5584SHong Zhang Level: intermediate 248d97c5584SHong Zhang 249d97c5584SHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation() 250d97c5584SHong Zhang @*/ 2517087cfbeSBarry Smith PetscErrorCode PetscSortSplit(PetscInt ncut,PetscInt n,PetscScalar a[],PetscInt idx[]) 252d97c5584SHong Zhang { 253d97c5584SHong Zhang PetscInt i,mid,last,itmp,j,first; 254d97c5584SHong Zhang PetscScalar d,tmp; 255d97c5584SHong Zhang PetscReal abskey; 256d97c5584SHong Zhang 257d97c5584SHong Zhang PetscFunctionBegin; 258d97c5584SHong Zhang first = 0; 259d97c5584SHong Zhang last = n-1; 260d97c5584SHong Zhang if (ncut < first || ncut > last) PetscFunctionReturn(0); 261d97c5584SHong Zhang 262d97c5584SHong Zhang while (1) { 263d97c5584SHong Zhang mid = first; 2642a274a78SSatish Balay d = a[mid]; 2652a274a78SSatish Balay abskey = PetscAbsScalar(d); 266d97c5584SHong Zhang i = last; 267d97c5584SHong Zhang for (j = first + 1; j <= i; ++j) { 2682a274a78SSatish Balay d = a[j]; 2692a274a78SSatish Balay if (PetscAbsScalar(d) >= abskey) { 270d97c5584SHong Zhang ++mid; 271d97c5584SHong Zhang /* interchange */ 272d97c5584SHong Zhang tmp = a[mid]; itmp = idx[mid]; 273d97c5584SHong Zhang a[mid] = a[j]; idx[mid] = idx[j]; 274d97c5584SHong Zhang a[j] = tmp; idx[j] = itmp; 275d97c5584SHong Zhang } 276d97c5584SHong Zhang } 277d97c5584SHong Zhang 278d97c5584SHong Zhang /* interchange */ 279d97c5584SHong Zhang tmp = a[mid]; itmp = idx[mid]; 280d97c5584SHong Zhang a[mid] = a[first]; idx[mid] = idx[first]; 281d97c5584SHong Zhang a[first] = tmp; idx[first] = itmp; 282d97c5584SHong Zhang 283d97c5584SHong Zhang /* test for while loop */ 284a297a907SKarl Rupp if (mid == ncut) break; 285a297a907SKarl Rupp else if (mid > ncut) last = mid - 1; 286a297a907SKarl Rupp else first = mid + 1; 287d97c5584SHong Zhang } 288d97c5584SHong Zhang PetscFunctionReturn(0); 289d97c5584SHong Zhang } 290021822bcSHong Zhang 291021822bcSHong Zhang /*@ 292021822bcSHong Zhang PetscSortSplitReal - Quick-sort split of an array of PetscReals in place. 293021822bcSHong Zhang 294021822bcSHong Zhang Not Collective 295021822bcSHong Zhang 296021822bcSHong Zhang Input Parameters: 297021822bcSHong Zhang + ncut - splitig index 298021822bcSHong Zhang . n - number of values to sort 299021822bcSHong Zhang . a - array of values in PetscReal 300021822bcSHong Zhang - idx - index for array a 301021822bcSHong Zhang 302021822bcSHong Zhang Output Parameters: 303021822bcSHong Zhang + a - permuted array of real values such that its elements satisfy: 304021822bcSHong Zhang abs(a[i]) >= abs(a[ncut-1]) for i < ncut and 305021822bcSHong Zhang abs(a[i]) <= abs(a[ncut-1]) for i >= ncut 306021822bcSHong Zhang - idx - permuted index of array a 307021822bcSHong Zhang 308021822bcSHong Zhang Level: intermediate 309021822bcSHong Zhang 310021822bcSHong Zhang .seealso: PetscSortInt(), PetscSortRealWithPermutation() 311021822bcSHong Zhang @*/ 3127087cfbeSBarry Smith PetscErrorCode PetscSortSplitReal(PetscInt ncut,PetscInt n,PetscReal a[],PetscInt idx[]) 313021822bcSHong Zhang { 314021822bcSHong Zhang PetscInt i,mid,last,itmp,j,first; 315021822bcSHong Zhang PetscReal d,tmp; 316021822bcSHong Zhang PetscReal abskey; 317021822bcSHong Zhang 318021822bcSHong Zhang PetscFunctionBegin; 319021822bcSHong Zhang first = 0; 320021822bcSHong Zhang last = n-1; 321021822bcSHong Zhang if (ncut < first || ncut > last) PetscFunctionReturn(0); 322021822bcSHong Zhang 323021822bcSHong Zhang while (1) { 324021822bcSHong Zhang mid = first; 3252a274a78SSatish Balay d = a[mid]; 3262a274a78SSatish Balay abskey = PetscAbsReal(d); 327021822bcSHong Zhang i = last; 328021822bcSHong Zhang for (j = first + 1; j <= i; ++j) { 3292a274a78SSatish Balay d = a[j]; 3302a274a78SSatish Balay if (PetscAbsReal(d) >= abskey) { 331021822bcSHong Zhang ++mid; 332021822bcSHong Zhang /* interchange */ 333021822bcSHong Zhang tmp = a[mid]; itmp = idx[mid]; 334021822bcSHong Zhang a[mid] = a[j]; idx[mid] = idx[j]; 335021822bcSHong Zhang a[j] = tmp; idx[j] = itmp; 336021822bcSHong Zhang } 337021822bcSHong Zhang } 338021822bcSHong Zhang 339021822bcSHong Zhang /* interchange */ 340021822bcSHong Zhang tmp = a[mid]; itmp = idx[mid]; 341021822bcSHong Zhang a[mid] = a[first]; idx[mid] = idx[first]; 342021822bcSHong Zhang a[first] = tmp; idx[first] = itmp; 343021822bcSHong Zhang 344021822bcSHong Zhang /* test for while loop */ 345a297a907SKarl Rupp if (mid == ncut) break; 346a297a907SKarl Rupp else if (mid > ncut) last = mid - 1; 347a297a907SKarl Rupp else first = mid + 1; 348021822bcSHong Zhang } 349021822bcSHong Zhang PetscFunctionReturn(0); 350021822bcSHong Zhang } 351021822bcSHong Zhang 352