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 5a5b23f4aSJose E. Roman of overhead in calling the comparison routines. 6e5c89e4eSSatish Balay 7e5c89e4eSSatish Balay */ 8c6db04a5SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 939f41f7fSStefano Zampini #include <petsc/private/petscimpl.h> 10e5c89e4eSSatish Balay 11*9371c9d4SSatish Balay #define SWAP(a, b, t) \ 12*9371c9d4SSatish Balay { \ 13*9371c9d4SSatish Balay t = a; \ 14*9371c9d4SSatish Balay a = b; \ 15*9371c9d4SSatish Balay b = t; \ 16*9371c9d4SSatish Balay } 17e5c89e4eSSatish Balay 186a7c706bSVaclav Hapla /*@ 196a7c706bSVaclav Hapla PetscSortedReal - Determines whether the array is sorted. 206a7c706bSVaclav Hapla 216a7c706bSVaclav Hapla Not Collective 226a7c706bSVaclav Hapla 236a7c706bSVaclav Hapla Input Parameters: 246a7c706bSVaclav Hapla + n - number of values 256a7c706bSVaclav Hapla - X - array of integers 266a7c706bSVaclav Hapla 276a7c706bSVaclav Hapla Output Parameters: 286a7c706bSVaclav Hapla . sorted - flag whether the array is sorted 296a7c706bSVaclav Hapla 306a7c706bSVaclav Hapla Level: intermediate 316a7c706bSVaclav Hapla 32db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortedInt()`, `PetscSortedMPIInt()` 336a7c706bSVaclav Hapla @*/ 34*9371c9d4SSatish Balay PetscErrorCode PetscSortedReal(PetscInt n, const PetscReal X[], PetscBool *sorted) { 356a7c706bSVaclav Hapla PetscFunctionBegin; 366a7c706bSVaclav Hapla PetscSorted(n, X, *sorted); 376a7c706bSVaclav Hapla PetscFunctionReturn(0); 386a7c706bSVaclav Hapla } 396a7c706bSVaclav Hapla 40e5c89e4eSSatish Balay /* A simple version of quicksort; taken from Kernighan and Ritchie, page 87 */ 41*9371c9d4SSatish Balay static PetscErrorCode PetscSortReal_Private(PetscReal *v, PetscInt right) { 42e5c89e4eSSatish Balay PetscInt i, last; 43e5c89e4eSSatish Balay PetscReal vl, tmp; 44e5c89e4eSSatish Balay 45e5c89e4eSSatish Balay PetscFunctionBegin; 46e5c89e4eSSatish Balay if (right <= 1) { 47e5c89e4eSSatish Balay if (right == 1) { 48e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP(v[0], v[1], tmp); 49e5c89e4eSSatish Balay } 50e5c89e4eSSatish Balay PetscFunctionReturn(0); 51e5c89e4eSSatish Balay } 52e5c89e4eSSatish Balay SWAP(v[0], v[right / 2], tmp); 53e5c89e4eSSatish Balay vl = v[0]; 54e5c89e4eSSatish Balay last = 0; 55e5c89e4eSSatish Balay for (i = 1; i <= right; i++) { 56*9371c9d4SSatish Balay if (v[i] < vl) { 57*9371c9d4SSatish Balay last++; 58*9371c9d4SSatish Balay SWAP(v[last], v[i], tmp); 59*9371c9d4SSatish Balay } 60e5c89e4eSSatish Balay } 61e5c89e4eSSatish Balay SWAP(v[0], v[last], tmp); 62e5c89e4eSSatish Balay PetscSortReal_Private(v, last - 1); 63e5c89e4eSSatish Balay PetscSortReal_Private(v + last + 1, right - (last + 1)); 64e5c89e4eSSatish Balay PetscFunctionReturn(0); 65e5c89e4eSSatish Balay } 66e5c89e4eSSatish Balay 67e5c89e4eSSatish Balay /*@ 68e5c89e4eSSatish Balay PetscSortReal - Sorts an array of doubles in place in increasing order. 69e5c89e4eSSatish Balay 70e5c89e4eSSatish Balay Not Collective 71e5c89e4eSSatish Balay 72e5c89e4eSSatish Balay Input Parameters: 73e5c89e4eSSatish Balay + n - number of values 74e5c89e4eSSatish Balay - v - array of doubles 75e5c89e4eSSatish Balay 76676f2a66SJacob Faibussowitsch Notes: 77676f2a66SJacob Faibussowitsch This function serves as an alternative to PetscRealSortSemiOrdered(), and may perform faster especially if the array 78a5b23f4aSJose E. Roman is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their 79676f2a66SJacob Faibussowitsch code to see which routine is fastest. 80676f2a66SJacob Faibussowitsch 81e5c89e4eSSatish Balay Level: intermediate 82e5c89e4eSSatish Balay 83db781477SPatrick Sanan .seealso: `PetscRealSortSemiOrdered()`, `PetscSortInt()`, `PetscSortRealWithPermutation()`, `PetscSortRealWithArrayInt()` 84e5c89e4eSSatish Balay @*/ 85*9371c9d4SSatish Balay PetscErrorCode PetscSortReal(PetscInt n, PetscReal v[]) { 86e5c89e4eSSatish Balay PetscInt j, k; 87e5c89e4eSSatish Balay PetscReal tmp, vk; 88e5c89e4eSSatish Balay 89e5c89e4eSSatish Balay PetscFunctionBegin; 90dadcf809SJacob Faibussowitsch PetscValidRealPointer(v, 2); 91e5c89e4eSSatish Balay if (n < 8) { 92e5c89e4eSSatish Balay for (k = 0; k < n; k++) { 93e5c89e4eSSatish Balay vk = v[k]; 94e5c89e4eSSatish Balay for (j = k + 1; j < n; j++) { 95e5c89e4eSSatish Balay if (vk > v[j]) { 96e5c89e4eSSatish Balay SWAP(v[k], v[j], tmp); 97e5c89e4eSSatish Balay vk = v[k]; 98e5c89e4eSSatish Balay } 99e5c89e4eSSatish Balay } 100e5c89e4eSSatish Balay } 101a297a907SKarl Rupp } else PetscSortReal_Private(v, n - 1); 102e5c89e4eSSatish Balay PetscFunctionReturn(0); 103e5c89e4eSSatish Balay } 104e5c89e4eSSatish Balay 105*9371c9d4SSatish Balay #define SWAP2ri(a, b, c, d, rt, it) \ 106*9371c9d4SSatish Balay { \ 107*9371c9d4SSatish Balay rt = a; \ 108*9371c9d4SSatish Balay a = b; \ 109*9371c9d4SSatish Balay b = rt; \ 110*9371c9d4SSatish Balay it = c; \ 111*9371c9d4SSatish Balay c = d; \ 112*9371c9d4SSatish Balay d = it; \ 113*9371c9d4SSatish Balay } 11439f41f7fSStefano Zampini 11539f41f7fSStefano Zampini /* modified from PetscSortIntWithArray_Private */ 116*9371c9d4SSatish Balay static PetscErrorCode PetscSortRealWithArrayInt_Private(PetscReal *v, PetscInt *V, PetscInt right) { 11739f41f7fSStefano Zampini PetscInt i, last, itmp; 11839f41f7fSStefano Zampini PetscReal rvl, rtmp; 11939f41f7fSStefano Zampini 12039f41f7fSStefano Zampini PetscFunctionBegin; 12139f41f7fSStefano Zampini if (right <= 1) { 12239f41f7fSStefano Zampini if (right == 1) { 12339f41f7fSStefano Zampini if (v[0] > v[1]) SWAP2ri(v[0], v[1], V[0], V[1], rtmp, itmp); 12439f41f7fSStefano Zampini } 12539f41f7fSStefano Zampini PetscFunctionReturn(0); 12639f41f7fSStefano Zampini } 12739f41f7fSStefano Zampini SWAP2ri(v[0], v[right / 2], V[0], V[right / 2], rtmp, itmp); 12839f41f7fSStefano Zampini rvl = v[0]; 12939f41f7fSStefano Zampini last = 0; 13039f41f7fSStefano Zampini for (i = 1; i <= right; i++) { 131*9371c9d4SSatish Balay if (v[i] < rvl) { 132*9371c9d4SSatish Balay last++; 133*9371c9d4SSatish Balay SWAP2ri(v[last], v[i], V[last], V[i], rtmp, itmp); 134*9371c9d4SSatish Balay } 13539f41f7fSStefano Zampini } 13639f41f7fSStefano Zampini SWAP2ri(v[0], v[last], V[0], V[last], rtmp, itmp); 1379566063dSJacob Faibussowitsch PetscCall(PetscSortRealWithArrayInt_Private(v, V, last - 1)); 1389566063dSJacob Faibussowitsch PetscCall(PetscSortRealWithArrayInt_Private(v + last + 1, V + last + 1, right - (last + 1))); 13939f41f7fSStefano Zampini PetscFunctionReturn(0); 14039f41f7fSStefano Zampini } 14139f41f7fSStefano Zampini /*@ 14239f41f7fSStefano Zampini PetscSortRealWithArrayInt - Sorts an array of PetscReal in place in increasing order; 14339f41f7fSStefano Zampini changes a second PetscInt array to match the sorted first array. 14439f41f7fSStefano Zampini 14539f41f7fSStefano Zampini Not Collective 14639f41f7fSStefano Zampini 14739f41f7fSStefano Zampini Input Parameters: 14839f41f7fSStefano Zampini + n - number of values 14939f41f7fSStefano Zampini . i - array of integers 15039f41f7fSStefano Zampini - I - second array of integers 15139f41f7fSStefano Zampini 15239f41f7fSStefano Zampini Level: intermediate 15339f41f7fSStefano Zampini 154db781477SPatrick Sanan .seealso: `PetscSortReal()` 15539f41f7fSStefano Zampini @*/ 156*9371c9d4SSatish Balay PetscErrorCode PetscSortRealWithArrayInt(PetscInt n, PetscReal r[], PetscInt Ii[]) { 15739f41f7fSStefano Zampini PetscInt j, k, itmp; 15839f41f7fSStefano Zampini PetscReal rk, rtmp; 15939f41f7fSStefano Zampini 16039f41f7fSStefano Zampini PetscFunctionBegin; 161dadcf809SJacob Faibussowitsch PetscValidRealPointer(r, 2); 162dadcf809SJacob Faibussowitsch PetscValidIntPointer(Ii, 3); 16339f41f7fSStefano Zampini if (n < 8) { 16439f41f7fSStefano Zampini for (k = 0; k < n; k++) { 16539f41f7fSStefano Zampini rk = r[k]; 16639f41f7fSStefano Zampini for (j = k + 1; j < n; j++) { 16739f41f7fSStefano Zampini if (rk > r[j]) { 16839f41f7fSStefano Zampini SWAP2ri(r[k], r[j], Ii[k], Ii[j], rtmp, itmp); 16939f41f7fSStefano Zampini rk = r[k]; 17039f41f7fSStefano Zampini } 17139f41f7fSStefano Zampini } 17239f41f7fSStefano Zampini } 17339f41f7fSStefano Zampini } else { 1749566063dSJacob Faibussowitsch PetscCall(PetscSortRealWithArrayInt_Private(r, Ii, n - 1)); 17539f41f7fSStefano Zampini } 17639f41f7fSStefano Zampini PetscFunctionReturn(0); 17739f41f7fSStefano Zampini } 17839f41f7fSStefano Zampini 17939f41f7fSStefano Zampini /*@ 18039f41f7fSStefano Zampini PetscFindReal - Finds a PetscReal in a sorted array of PetscReals 18139f41f7fSStefano Zampini 18239f41f7fSStefano Zampini Not Collective 18339f41f7fSStefano Zampini 18439f41f7fSStefano Zampini Input Parameters: 18539f41f7fSStefano Zampini + key - the value to locate 18639f41f7fSStefano Zampini . n - number of values in the array 18739f41f7fSStefano Zampini . ii - array of values 18839f41f7fSStefano Zampini - eps - tolerance used to compare 18939f41f7fSStefano Zampini 19039f41f7fSStefano Zampini Output Parameter: 19139f41f7fSStefano Zampini . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 19239f41f7fSStefano Zampini 19339f41f7fSStefano Zampini Level: intermediate 19439f41f7fSStefano Zampini 195db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortRealWithArrayInt()` 19639f41f7fSStefano Zampini @*/ 197*9371c9d4SSatish Balay PetscErrorCode PetscFindReal(PetscReal key, PetscInt n, const PetscReal t[], PetscReal eps, PetscInt *loc) { 19839f41f7fSStefano Zampini PetscInt lo = 0, hi = n; 19939f41f7fSStefano Zampini 20039f41f7fSStefano Zampini PetscFunctionBegin; 201dadcf809SJacob Faibussowitsch PetscValidIntPointer(loc, 5); 202*9371c9d4SSatish Balay if (!n) { 203*9371c9d4SSatish Balay *loc = -1; 204*9371c9d4SSatish Balay PetscFunctionReturn(0); 205*9371c9d4SSatish Balay } 206dadcf809SJacob Faibussowitsch PetscValidRealPointer(t, 3); 2076a7c706bSVaclav Hapla PetscCheckSorted(n, t); 20839f41f7fSStefano Zampini while (hi - lo > 1) { 20939f41f7fSStefano Zampini PetscInt mid = lo + (hi - lo) / 2; 21039f41f7fSStefano Zampini if (key < t[mid]) hi = mid; 21139f41f7fSStefano Zampini else lo = mid; 21239f41f7fSStefano Zampini } 21339f41f7fSStefano Zampini *loc = (PetscAbsReal(key - t[lo]) < eps) ? lo : -(lo + (key > t[lo]) + 1); 21439f41f7fSStefano Zampini PetscFunctionReturn(0); 21539f41f7fSStefano Zampini } 216745b41b2SMatthew G. Knepley 217745b41b2SMatthew G. Knepley /*@ 218745b41b2SMatthew G. Knepley PetscSortRemoveDupsReal - Sorts an array of doubles in place in increasing order removes all duplicate entries 219745b41b2SMatthew G. Knepley 220745b41b2SMatthew G. Knepley Not Collective 221745b41b2SMatthew G. Knepley 222745b41b2SMatthew G. Knepley Input Parameters: 223745b41b2SMatthew G. Knepley + n - number of values 224745b41b2SMatthew G. Knepley - v - array of doubles 225745b41b2SMatthew G. Knepley 226745b41b2SMatthew G. Knepley Output Parameter: 227745b41b2SMatthew G. Knepley . n - number of non-redundant values 228745b41b2SMatthew G. Knepley 229745b41b2SMatthew G. Knepley Level: intermediate 230745b41b2SMatthew G. Knepley 231db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortRemoveDupsInt()` 232745b41b2SMatthew G. Knepley @*/ 233*9371c9d4SSatish Balay PetscErrorCode PetscSortRemoveDupsReal(PetscInt *n, PetscReal v[]) { 234745b41b2SMatthew G. Knepley PetscInt i, s = 0, N = *n, b = 0; 235745b41b2SMatthew G. Knepley 236745b41b2SMatthew G. Knepley PetscFunctionBegin; 2379566063dSJacob Faibussowitsch PetscCall(PetscSortReal(N, v)); 238745b41b2SMatthew G. Knepley for (i = 0; i < N - 1; i++) { 239745b41b2SMatthew G. Knepley if (v[b + s + 1] != v[b]) { 240*9371c9d4SSatish Balay v[b + 1] = v[b + s + 1]; 241*9371c9d4SSatish Balay b++; 242745b41b2SMatthew G. Knepley } else s++; 243745b41b2SMatthew G. Knepley } 244745b41b2SMatthew G. Knepley *n = N - s; 245745b41b2SMatthew G. Knepley PetscFunctionReturn(0); 246745b41b2SMatthew G. Knepley } 247745b41b2SMatthew G. Knepley 248d97c5584SHong Zhang /*@ 249021822bcSHong Zhang PetscSortSplit - Quick-sort split of an array of PetscScalars in place. 250d97c5584SHong Zhang 251d97c5584SHong Zhang Not Collective 252d97c5584SHong Zhang 253d97c5584SHong Zhang Input Parameters: 254d97c5584SHong Zhang + ncut - splitig index 2556b867d5aSJose E. Roman - n - number of values to sort 256d97c5584SHong Zhang 2576b867d5aSJose E. Roman Input/Output Parameters: 2586b867d5aSJose E. Roman + a - array of values, on output the values are permuted such that its elements satisfy: 259d97c5584SHong Zhang abs(a[i]) >= abs(a[ncut-1]) for i < ncut and 260d97c5584SHong Zhang abs(a[i]) <= abs(a[ncut-1]) for i >= ncut 2616b867d5aSJose E. Roman - idx - index for array a, on output permuted accordingly 262d97c5584SHong Zhang 263d97c5584SHong Zhang Level: intermediate 264d97c5584SHong Zhang 265db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()` 266d97c5584SHong Zhang @*/ 267*9371c9d4SSatish Balay PetscErrorCode PetscSortSplit(PetscInt ncut, PetscInt n, PetscScalar a[], PetscInt idx[]) { 268d97c5584SHong Zhang PetscInt i, mid, last, itmp, j, first; 269d97c5584SHong Zhang PetscScalar d, tmp; 270d97c5584SHong Zhang PetscReal abskey; 271d97c5584SHong Zhang 272d97c5584SHong Zhang PetscFunctionBegin; 273d97c5584SHong Zhang first = 0; 274d97c5584SHong Zhang last = n - 1; 275d97c5584SHong Zhang if (ncut < first || ncut > last) PetscFunctionReturn(0); 276d97c5584SHong Zhang 277d97c5584SHong Zhang while (1) { 278d97c5584SHong Zhang mid = first; 2792a274a78SSatish Balay d = a[mid]; 2802a274a78SSatish Balay abskey = PetscAbsScalar(d); 281d97c5584SHong Zhang i = last; 282d97c5584SHong Zhang for (j = first + 1; j <= i; ++j) { 2832a274a78SSatish Balay d = a[j]; 2842a274a78SSatish Balay if (PetscAbsScalar(d) >= abskey) { 285d97c5584SHong Zhang ++mid; 286d97c5584SHong Zhang /* interchange */ 287*9371c9d4SSatish Balay tmp = a[mid]; 288*9371c9d4SSatish Balay itmp = idx[mid]; 289*9371c9d4SSatish Balay a[mid] = a[j]; 290*9371c9d4SSatish Balay idx[mid] = idx[j]; 291*9371c9d4SSatish Balay a[j] = tmp; 292*9371c9d4SSatish Balay idx[j] = itmp; 293d97c5584SHong Zhang } 294d97c5584SHong Zhang } 295d97c5584SHong Zhang 296d97c5584SHong Zhang /* interchange */ 297*9371c9d4SSatish Balay tmp = a[mid]; 298*9371c9d4SSatish Balay itmp = idx[mid]; 299*9371c9d4SSatish Balay a[mid] = a[first]; 300*9371c9d4SSatish Balay idx[mid] = idx[first]; 301*9371c9d4SSatish Balay a[first] = tmp; 302*9371c9d4SSatish Balay idx[first] = itmp; 303d97c5584SHong Zhang 304d97c5584SHong Zhang /* test for while loop */ 305a297a907SKarl Rupp if (mid == ncut) break; 306a297a907SKarl Rupp else if (mid > ncut) last = mid - 1; 307a297a907SKarl Rupp else first = mid + 1; 308d97c5584SHong Zhang } 309d97c5584SHong Zhang PetscFunctionReturn(0); 310d97c5584SHong Zhang } 311021822bcSHong Zhang 312021822bcSHong Zhang /*@ 313021822bcSHong Zhang PetscSortSplitReal - Quick-sort split of an array of PetscReals in place. 314021822bcSHong Zhang 315021822bcSHong Zhang Not Collective 316021822bcSHong Zhang 317021822bcSHong Zhang Input Parameters: 318021822bcSHong Zhang + ncut - splitig index 3196b867d5aSJose E. Roman - n - number of values to sort 320021822bcSHong Zhang 3216b867d5aSJose E. Roman Input/Output Parameters: 3226b867d5aSJose E. Roman + a - array of values, on output the values are permuted such that its elements satisfy: 323021822bcSHong Zhang abs(a[i]) >= abs(a[ncut-1]) for i < ncut and 324021822bcSHong Zhang abs(a[i]) <= abs(a[ncut-1]) for i >= ncut 3256b867d5aSJose E. Roman - idx - index for array a, on output permuted accordingly 326021822bcSHong Zhang 327021822bcSHong Zhang Level: intermediate 328021822bcSHong Zhang 329db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()` 330021822bcSHong Zhang @*/ 331*9371c9d4SSatish Balay PetscErrorCode PetscSortSplitReal(PetscInt ncut, PetscInt n, PetscReal a[], PetscInt idx[]) { 332021822bcSHong Zhang PetscInt i, mid, last, itmp, j, first; 333021822bcSHong Zhang PetscReal d, tmp; 334021822bcSHong Zhang PetscReal abskey; 335021822bcSHong Zhang 336021822bcSHong Zhang PetscFunctionBegin; 337021822bcSHong Zhang first = 0; 338021822bcSHong Zhang last = n - 1; 339021822bcSHong Zhang if (ncut < first || ncut > last) PetscFunctionReturn(0); 340021822bcSHong Zhang 341021822bcSHong Zhang while (1) { 342021822bcSHong Zhang mid = first; 3432a274a78SSatish Balay d = a[mid]; 3442a274a78SSatish Balay abskey = PetscAbsReal(d); 345021822bcSHong Zhang i = last; 346021822bcSHong Zhang for (j = first + 1; j <= i; ++j) { 3472a274a78SSatish Balay d = a[j]; 3482a274a78SSatish Balay if (PetscAbsReal(d) >= abskey) { 349021822bcSHong Zhang ++mid; 350021822bcSHong Zhang /* interchange */ 351*9371c9d4SSatish Balay tmp = a[mid]; 352*9371c9d4SSatish Balay itmp = idx[mid]; 353*9371c9d4SSatish Balay a[mid] = a[j]; 354*9371c9d4SSatish Balay idx[mid] = idx[j]; 355*9371c9d4SSatish Balay a[j] = tmp; 356*9371c9d4SSatish Balay idx[j] = itmp; 357021822bcSHong Zhang } 358021822bcSHong Zhang } 359021822bcSHong Zhang 360021822bcSHong Zhang /* interchange */ 361*9371c9d4SSatish Balay tmp = a[mid]; 362*9371c9d4SSatish Balay itmp = idx[mid]; 363*9371c9d4SSatish Balay a[mid] = a[first]; 364*9371c9d4SSatish Balay idx[mid] = idx[first]; 365*9371c9d4SSatish Balay a[first] = tmp; 366*9371c9d4SSatish Balay idx[first] = itmp; 367021822bcSHong Zhang 368021822bcSHong Zhang /* test for while loop */ 369a297a907SKarl Rupp if (mid == ncut) break; 370a297a907SKarl Rupp else if (mid > ncut) last = mid - 1; 371a297a907SKarl Rupp else first = mid + 1; 372021822bcSHong Zhang } 373021822bcSHong Zhang PetscFunctionReturn(0); 374021822bcSHong Zhang } 375