17d0a6c19SBarry Smith 2e5c89e4eSSatish Balay /* 3e5c89e4eSSatish Balay This file contains routines for sorting integers. Values are sorted in place. 4c4762a1bSJed Brown One can use src/sys/tests/ex52.c for benchmarking. 5e5c89e4eSSatish Balay */ 6af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I "petscsys.h" I*/ 7f1cab4e1SJunchao Zhang #include <petsc/private/hashseti.h> 8e5c89e4eSSatish Balay 9a8582168SJed Brown #define MEDIAN3(v,a,b,c) \ 10a8582168SJed Brown (v[a]<v[b] \ 11a8582168SJed Brown ? (v[b]<v[c] \ 1259e16d97SJunchao Zhang ? (b) \ 1359e16d97SJunchao Zhang : (v[a]<v[c] ? (c) : (a))) \ 14a8582168SJed Brown : (v[c]<v[b] \ 1559e16d97SJunchao Zhang ? (b) \ 1659e16d97SJunchao Zhang : (v[a]<v[c] ? (a) : (c)))) 17a8582168SJed Brown 18a8582168SJed Brown #define MEDIAN(v,right) MEDIAN3(v,right/4,right/2,right/4*3) 19ef8e3583SJed Brown 2059e16d97SJunchao Zhang /* Swap one, two or three pairs. Each pair can have its own type */ 2159e16d97SJunchao Zhang #define SWAP1(a,b,t1) do {t1=a;a=b;b=t1;} while (0) 2259e16d97SJunchao Zhang #define SWAP2(a,b,c,d,t1,t2) do {t1=a;a=b;b=t1; t2=c;c=d;d=t2;} while (0) 2359e16d97SJunchao Zhang #define SWAP3(a,b,c,d,e,f,t1,t2,t3) do {t1=a;a=b;b=t1; t2=c;c=d;d=t2; t3=e;e=f;f=t3;} while (0) 2459e16d97SJunchao Zhang 2559e16d97SJunchao Zhang /* Swap a & b, *c & *d. c, d, t2 are pointers to a type of size <siz> */ 2659e16d97SJunchao Zhang #define SWAP2Data(a,b,c,d,t1,t2,siz) \ 2759e16d97SJunchao Zhang do { \ 2859e16d97SJunchao Zhang PetscErrorCode ierr; \ 2959e16d97SJunchao Zhang t1=a; a=b; b=t1; \ 3059e16d97SJunchao Zhang ierr = PetscMemcpy(t2,c,siz);CHKERRQ(ierr); \ 3159e16d97SJunchao Zhang ierr = PetscMemcpy(c,d,siz);CHKERRQ(ierr); \ 3259e16d97SJunchao Zhang ierr = PetscMemcpy(d,t2,siz);CHKERRQ(ierr); \ 3359e16d97SJunchao Zhang } while (0) 34e5c89e4eSSatish Balay 35e5c89e4eSSatish Balay /* 3659e16d97SJunchao Zhang Partition X[lo,hi] into two parts: X[lo,l) <= pivot; X[r,hi] > pivot 37e5c89e4eSSatish Balay 3859e16d97SJunchao Zhang Input Parameters: 3959e16d97SJunchao Zhang + X - array to partition 4059e16d97SJunchao Zhang . pivot - a pivot of X[] 4159e16d97SJunchao Zhang . t1 - temp variable for X 4259e16d97SJunchao Zhang - lo,hi - lower and upper bound of the array 4359e16d97SJunchao Zhang 4459e16d97SJunchao Zhang Output Parameters: 4559e16d97SJunchao Zhang + l,r - of type PetscInt 4659e16d97SJunchao Zhang 4759e16d97SJunchao Zhang Notes: 4859e16d97SJunchao Zhang The TwoWayPartition2/3 variants also partition other arrays along with X. 4959e16d97SJunchao Zhang These arrays can have different types, so they provide their own temp t2,t3 5059e16d97SJunchao Zhang */ 5159e16d97SJunchao Zhang #define TwoWayPartition1(X,pivot,t1,lo,hi,l,r) \ 5259e16d97SJunchao Zhang do { \ 5359e16d97SJunchao Zhang l = lo; \ 5459e16d97SJunchao Zhang r = hi; \ 5559e16d97SJunchao Zhang while (1) { \ 5659e16d97SJunchao Zhang while (X[l] < pivot) l++; \ 5759e16d97SJunchao Zhang while (X[r] > pivot) r--; \ 5859e16d97SJunchao Zhang if (l >= r) {r++; break;} \ 5959e16d97SJunchao Zhang SWAP1(X[l],X[r],t1); \ 6059e16d97SJunchao Zhang l++; \ 6159e16d97SJunchao Zhang r--; \ 6259e16d97SJunchao Zhang } \ 6359e16d97SJunchao Zhang } while (0) 6459e16d97SJunchao Zhang 65ce605777SToby Isaac /* 66ce605777SToby Isaac Partition X[lo,hi] into two parts: X[lo,l) >= pivot; X[r,hi] < pivot 67ce605777SToby Isaac 68ce605777SToby Isaac Input Parameters: 69ce605777SToby Isaac + X - array to partition 70ce605777SToby Isaac . pivot - a pivot of X[] 71ce605777SToby Isaac . t1 - temp variable for X 72ce605777SToby Isaac - lo,hi - lower and upper bound of the array 73ce605777SToby Isaac 74ce605777SToby Isaac Output Parameters: 75ce605777SToby Isaac + l,r - of type PetscInt 76ce605777SToby Isaac 77ce605777SToby Isaac Notes: 78ce605777SToby Isaac The TwoWayPartition2/3 variants also partition other arrays along with X. 79ce605777SToby Isaac These arrays can have different types, so they provide their own temp t2,t3 80ce605777SToby Isaac */ 81ce605777SToby Isaac #define TwoWayPartitionReverse1(X,pivot,t1,lo,hi,l,r) \ 82ce605777SToby Isaac do { \ 83ce605777SToby Isaac l = lo; \ 84ce605777SToby Isaac r = hi; \ 85ce605777SToby Isaac while (1) { \ 86ce605777SToby Isaac while (X[l] > pivot) l++; \ 87ce605777SToby Isaac while (X[r] < pivot) r--; \ 88ce605777SToby Isaac if (l >= r) {r++; break;} \ 89ce605777SToby Isaac SWAP1(X[l],X[r],t1); \ 90ce605777SToby Isaac l++; \ 91ce605777SToby Isaac r--; \ 92ce605777SToby Isaac } \ 93ce605777SToby Isaac } while (0) 94ce605777SToby Isaac 9559e16d97SJunchao Zhang #define TwoWayPartition2(X,Y,pivot,t1,t2,lo,hi,l,r) \ 9659e16d97SJunchao Zhang do { \ 9759e16d97SJunchao Zhang l = lo; \ 9859e16d97SJunchao Zhang r = hi; \ 9959e16d97SJunchao Zhang while (1) { \ 10059e16d97SJunchao Zhang while (X[l] < pivot) l++; \ 10159e16d97SJunchao Zhang while (X[r] > pivot) r--; \ 10259e16d97SJunchao Zhang if (l >= r) {r++; break;} \ 10359e16d97SJunchao Zhang SWAP2(X[l],X[r],Y[l],Y[r],t1,t2); \ 10459e16d97SJunchao Zhang l++; \ 10559e16d97SJunchao Zhang r--; \ 10659e16d97SJunchao Zhang } \ 10759e16d97SJunchao Zhang } while (0) 10859e16d97SJunchao Zhang 10959e16d97SJunchao Zhang #define TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,lo,hi,l,r) \ 11059e16d97SJunchao Zhang do { \ 11159e16d97SJunchao Zhang l = lo; \ 11259e16d97SJunchao Zhang r = hi; \ 11359e16d97SJunchao Zhang while (1) { \ 11459e16d97SJunchao Zhang while (X[l] < pivot) l++; \ 11559e16d97SJunchao Zhang while (X[r] > pivot) r--; \ 11659e16d97SJunchao Zhang if (l >= r) {r++; break;} \ 11759e16d97SJunchao Zhang SWAP3(X[l],X[r],Y[l],Y[r],Z[l],Z[r],t1,t2,t3); \ 11859e16d97SJunchao Zhang l++; \ 11959e16d97SJunchao Zhang r--; \ 12059e16d97SJunchao Zhang } \ 12159e16d97SJunchao Zhang } while (0) 12259e16d97SJunchao Zhang 12359e16d97SJunchao Zhang /* Templates for similar functions used below */ 12459e16d97SJunchao Zhang #define QuickSort1(FuncName,X,n,pivot,t1,ierr) \ 12559e16d97SJunchao Zhang do { \ 12659e16d97SJunchao Zhang PetscInt i,j,p,l,r,hi=n-1; \ 12759e16d97SJunchao Zhang if (n < 8) { \ 12859e16d97SJunchao Zhang for (i=0; i<n; i++) { \ 12959e16d97SJunchao Zhang pivot = X[i]; \ 13059e16d97SJunchao Zhang for (j=i+1; j<n; j++) { \ 13159e16d97SJunchao Zhang if (pivot > X[j]) { \ 13259e16d97SJunchao Zhang SWAP1(X[i],X[j],t1); \ 13359e16d97SJunchao Zhang pivot = X[i]; \ 13459e16d97SJunchao Zhang } \ 13559e16d97SJunchao Zhang } \ 13659e16d97SJunchao Zhang } \ 13759e16d97SJunchao Zhang } else { \ 13859e16d97SJunchao Zhang p = MEDIAN(X,hi); \ 13959e16d97SJunchao Zhang pivot = X[p]; \ 14059e16d97SJunchao Zhang TwoWayPartition1(X,pivot,t1,0,hi,l,r); \ 14159e16d97SJunchao Zhang ierr = FuncName(l,X);CHKERRQ(ierr); \ 14259e16d97SJunchao Zhang ierr = FuncName(hi-r+1,X+r);CHKERRQ(ierr); \ 14359e16d97SJunchao Zhang } \ 14459e16d97SJunchao Zhang } while (0) 14559e16d97SJunchao Zhang 146ce605777SToby Isaac /* Templates for similar functions used below */ 147ce605777SToby Isaac #define QuickSortReverse1(FuncName,X,n,pivot,t1,ierr) \ 148ce605777SToby Isaac do { \ 149ce605777SToby Isaac PetscInt i,j,p,l,r,hi=n-1; \ 150ce605777SToby Isaac if (n < 8) { \ 151ce605777SToby Isaac for (i=0; i<n; i++) { \ 152ce605777SToby Isaac pivot = X[i]; \ 153ce605777SToby Isaac for (j=i+1; j<n; j++) { \ 154ce605777SToby Isaac if (pivot < X[j]) { \ 155ce605777SToby Isaac SWAP1(X[i],X[j],t1); \ 156ce605777SToby Isaac pivot = X[i]; \ 157ce605777SToby Isaac } \ 158ce605777SToby Isaac } \ 159ce605777SToby Isaac } \ 160ce605777SToby Isaac } else { \ 161ce605777SToby Isaac p = MEDIAN(X,hi); \ 162ce605777SToby Isaac pivot = X[p]; \ 163ce605777SToby Isaac TwoWayPartitionReverse1(X,pivot,t1,0,hi,l,r); \ 164ce605777SToby Isaac ierr = FuncName(l,X);CHKERRQ(ierr); \ 165ce605777SToby Isaac ierr = FuncName(hi-r+1,X+r);CHKERRQ(ierr); \ 166ce605777SToby Isaac } \ 167ce605777SToby Isaac } while (0) 168ce605777SToby Isaac 16959e16d97SJunchao Zhang #define QuickSort2(FuncName,X,Y,n,pivot,t1,t2,ierr) \ 17059e16d97SJunchao Zhang do { \ 17159e16d97SJunchao Zhang PetscInt i,j,p,l,r,hi=n-1; \ 17259e16d97SJunchao Zhang if (n < 8) { \ 17359e16d97SJunchao Zhang for (i=0; i<n; i++) { \ 17459e16d97SJunchao Zhang pivot = X[i]; \ 17559e16d97SJunchao Zhang for (j=i+1; j<n; j++) { \ 17659e16d97SJunchao Zhang if (pivot > X[j]) { \ 17759e16d97SJunchao Zhang SWAP2(X[i],X[j],Y[i],Y[j],t1,t2); \ 17859e16d97SJunchao Zhang pivot = X[i]; \ 17959e16d97SJunchao Zhang } \ 18059e16d97SJunchao Zhang } \ 18159e16d97SJunchao Zhang } \ 18259e16d97SJunchao Zhang } else { \ 18359e16d97SJunchao Zhang p = MEDIAN(X,hi); \ 18459e16d97SJunchao Zhang pivot = X[p]; \ 18559e16d97SJunchao Zhang TwoWayPartition2(X,Y,pivot,t1,t2,0,hi,l,r); \ 18659e16d97SJunchao Zhang ierr = FuncName(l,X,Y);CHKERRQ(ierr); \ 18759e16d97SJunchao Zhang ierr = FuncName(hi-r+1,X+r,Y+r);CHKERRQ(ierr); \ 18859e16d97SJunchao Zhang } \ 18959e16d97SJunchao Zhang } while (0) 19059e16d97SJunchao Zhang 19159e16d97SJunchao Zhang #define QuickSort3(FuncName,X,Y,Z,n,pivot,t1,t2,t3,ierr) \ 19259e16d97SJunchao Zhang do { \ 19359e16d97SJunchao Zhang PetscInt i,j,p,l,r,hi=n-1; \ 19459e16d97SJunchao Zhang if (n < 8) { \ 19559e16d97SJunchao Zhang for (i=0; i<n; i++) { \ 19659e16d97SJunchao Zhang pivot = X[i]; \ 19759e16d97SJunchao Zhang for (j=i+1; j<n; j++) { \ 19859e16d97SJunchao Zhang if (pivot > X[j]) { \ 19959e16d97SJunchao Zhang SWAP3(X[i],X[j],Y[i],Y[j],Z[i],Z[j],t1,t2,t3); \ 20059e16d97SJunchao Zhang pivot = X[i]; \ 20159e16d97SJunchao Zhang } \ 20259e16d97SJunchao Zhang } \ 20359e16d97SJunchao Zhang } \ 20459e16d97SJunchao Zhang } else { \ 20559e16d97SJunchao Zhang p = MEDIAN(X,hi); \ 20659e16d97SJunchao Zhang pivot = X[p]; \ 20759e16d97SJunchao Zhang TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,0,hi,l,r); \ 20859e16d97SJunchao Zhang ierr = FuncName(l,X,Y,Z);CHKERRQ(ierr); \ 20959e16d97SJunchao Zhang ierr = FuncName(hi-r+1,X+r,Y+r,Z+r);CHKERRQ(ierr); \ 21059e16d97SJunchao Zhang } \ 21159e16d97SJunchao Zhang } while (0) 212e5c89e4eSSatish Balay 213e5c89e4eSSatish Balay /*@ 2146a7c706bSVaclav Hapla PetscSortedInt - Determines whether the array is sorted. 2156a7c706bSVaclav Hapla 2166a7c706bSVaclav Hapla Not Collective 2176a7c706bSVaclav Hapla 2186a7c706bSVaclav Hapla Input Parameters: 2196a7c706bSVaclav Hapla + n - number of values 2206a7c706bSVaclav Hapla - X - array of integers 2216a7c706bSVaclav Hapla 2226a7c706bSVaclav Hapla Output Parameters: 2236a7c706bSVaclav Hapla . sorted - flag whether the array is sorted 2246a7c706bSVaclav Hapla 2256a7c706bSVaclav Hapla Level: intermediate 2266a7c706bSVaclav Hapla 2276a7c706bSVaclav Hapla .seealso: PetscSortInt(), PetscSortedMPIInt(), PetscSortedReal() 2286a7c706bSVaclav Hapla @*/ 2296a7c706bSVaclav Hapla PetscErrorCode PetscSortedInt(PetscInt n,const PetscInt X[],PetscBool *sorted) 2306a7c706bSVaclav Hapla { 2316a7c706bSVaclav Hapla PetscFunctionBegin; 2326a7c706bSVaclav Hapla PetscSorted(n,X,*sorted); 2336a7c706bSVaclav Hapla PetscFunctionReturn(0); 2346a7c706bSVaclav Hapla } 2356a7c706bSVaclav Hapla 2366a7c706bSVaclav Hapla /*@ 237e5c89e4eSSatish Balay PetscSortInt - Sorts an array of integers in place in increasing order. 238e5c89e4eSSatish Balay 239e5c89e4eSSatish Balay Not Collective 240e5c89e4eSSatish Balay 241e5c89e4eSSatish Balay Input Parameters: 242e5c89e4eSSatish Balay + n - number of values 24359e16d97SJunchao Zhang - X - array of integers 244e5c89e4eSSatish Balay 245*676f2a66SJacob Faibussowitsch Notes: 246*676f2a66SJacob Faibussowitsch This function serves as an alternative to PetscIntSortSemiOrdered(), and may perform faster especially if the array 247*676f2a66SJacob Faibussowitsch is completely random. There are exceptions to this and so it is __highly__ recomended that the user benchmark their 248*676f2a66SJacob Faibussowitsch code to see which routine is fastest. 249*676f2a66SJacob Faibussowitsch 250e5c89e4eSSatish Balay Level: intermediate 251e5c89e4eSSatish Balay 252*676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation() 253e5c89e4eSSatish Balay @*/ 254b4f531b9SJunchao Zhang PetscErrorCode PetscSortInt(PetscInt n,PetscInt X[]) 255e5c89e4eSSatish Balay { 25659e16d97SJunchao Zhang PetscErrorCode ierr; 25759e16d97SJunchao Zhang PetscInt pivot,t1; 258e5c89e4eSSatish Balay 259e5c89e4eSSatish Balay PetscFunctionBegin; 26059e16d97SJunchao Zhang QuickSort1(PetscSortInt,X,n,pivot,t1,ierr); 261e5c89e4eSSatish Balay PetscFunctionReturn(0); 262e5c89e4eSSatish Balay } 263e5c89e4eSSatish Balay 2641db36a52SBarry Smith /*@ 265ce605777SToby Isaac PetscSortReverseInt - Sorts an array of integers in place in decreasing order. 266ce605777SToby Isaac 267ce605777SToby Isaac Not Collective 268ce605777SToby Isaac 269ce605777SToby Isaac Input Parameters: 270ce605777SToby Isaac + n - number of values 271ce605777SToby Isaac - X - array of integers 272ce605777SToby Isaac 273ce605777SToby Isaac Level: intermediate 274ce605777SToby Isaac 275*676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithPermutation() 276ce605777SToby Isaac @*/ 277ce605777SToby Isaac PetscErrorCode PetscSortReverseInt(PetscInt n,PetscInt X[]) 278ce605777SToby Isaac { 279ce605777SToby Isaac PetscErrorCode ierr; 280ce605777SToby Isaac PetscInt pivot,t1; 281ce605777SToby Isaac 282ce605777SToby Isaac PetscFunctionBegin; 283ce605777SToby Isaac QuickSortReverse1(PetscSortReverseInt,X,n,pivot,t1,ierr); 284ce605777SToby Isaac PetscFunctionReturn(0); 285ce605777SToby Isaac } 286ce605777SToby Isaac 287ce605777SToby Isaac /*@ 28822ab5688SLisandro Dalcin PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted input array 28922ab5688SLisandro Dalcin 29022ab5688SLisandro Dalcin Not Collective 29122ab5688SLisandro Dalcin 29222ab5688SLisandro Dalcin Input Parameters: 29322ab5688SLisandro Dalcin + n - number of values 29459e16d97SJunchao Zhang - X - sorted array of integers 29522ab5688SLisandro Dalcin 29622ab5688SLisandro Dalcin Output Parameter: 29722ab5688SLisandro Dalcin . n - number of non-redundant values 29822ab5688SLisandro Dalcin 29922ab5688SLisandro Dalcin Level: intermediate 30022ab5688SLisandro Dalcin 30122ab5688SLisandro Dalcin .seealso: PetscSortInt() 30222ab5688SLisandro Dalcin @*/ 30359e16d97SJunchao Zhang PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n,PetscInt X[]) 30422ab5688SLisandro Dalcin { 30522ab5688SLisandro Dalcin PetscInt i,s = 0,N = *n, b = 0; 30622ab5688SLisandro Dalcin 30722ab5688SLisandro Dalcin PetscFunctionBegin; 308fb61b9e4SStefano Zampini PetscCheckSorted(*n,X); 30922ab5688SLisandro Dalcin for (i=0; i<N-1; i++) { 31059e16d97SJunchao Zhang if (X[b+s+1] != X[b]) { 31159e16d97SJunchao Zhang X[b+1] = X[b+s+1]; b++; 31222ab5688SLisandro Dalcin } else s++; 31322ab5688SLisandro Dalcin } 31422ab5688SLisandro Dalcin *n = N - s; 31522ab5688SLisandro Dalcin PetscFunctionReturn(0); 31622ab5688SLisandro Dalcin } 31722ab5688SLisandro Dalcin 31822ab5688SLisandro Dalcin /*@ 3191db36a52SBarry Smith PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries 3201db36a52SBarry Smith 3211db36a52SBarry Smith Not Collective 3221db36a52SBarry Smith 3231db36a52SBarry Smith Input Parameters: 3241db36a52SBarry Smith + n - number of values 32559e16d97SJunchao Zhang - X - array of integers 3261db36a52SBarry Smith 3271db36a52SBarry Smith Output Parameter: 3281db36a52SBarry Smith . n - number of non-redundant values 3291db36a52SBarry Smith 3301db36a52SBarry Smith Level: intermediate 3311db36a52SBarry Smith 332*676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortedRemoveDupsInt() 3331db36a52SBarry Smith @*/ 33459e16d97SJunchao Zhang PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n,PetscInt X[]) 3351db36a52SBarry Smith { 3361db36a52SBarry Smith PetscErrorCode ierr; 3371db36a52SBarry Smith 3381db36a52SBarry Smith PetscFunctionBegin; 33959e16d97SJunchao Zhang ierr = PetscSortInt(*n,X);CHKERRQ(ierr); 34059e16d97SJunchao Zhang ierr = PetscSortedRemoveDupsInt(n,X);CHKERRQ(ierr); 3411db36a52SBarry Smith PetscFunctionReturn(0); 3421db36a52SBarry Smith } 3431db36a52SBarry Smith 34460e03357SMatthew G Knepley /*@ 345d9f1162dSJed Brown PetscFindInt - Finds integer in a sorted array of integers 34660e03357SMatthew G Knepley 34760e03357SMatthew G Knepley Not Collective 34860e03357SMatthew G Knepley 34960e03357SMatthew G Knepley Input Parameters: 35060e03357SMatthew G Knepley + key - the integer to locate 351d9f1162dSJed Brown . n - number of values in the array 35259e16d97SJunchao Zhang - X - array of integers 35360e03357SMatthew G Knepley 35460e03357SMatthew G Knepley Output Parameter: 355d9f1162dSJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 35660e03357SMatthew G Knepley 35760e03357SMatthew G Knepley Level: intermediate 35860e03357SMatthew G Knepley 359*676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt() 36060e03357SMatthew G Knepley @*/ 36159e16d97SJunchao Zhang PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc) 36260e03357SMatthew G Knepley { 363d9f1162dSJed Brown PetscInt lo = 0,hi = n; 36460e03357SMatthew G Knepley 36560e03357SMatthew G Knepley PetscFunctionBegin; 36660e03357SMatthew G Knepley PetscValidPointer(loc,4); 36798781d82SMatthew G Knepley if (!n) {*loc = -1; PetscFunctionReturn(0);} 36859e16d97SJunchao Zhang PetscValidPointer(X,3); 3696a7c706bSVaclav Hapla PetscCheckSorted(n,X); 370d9f1162dSJed Brown while (hi - lo > 1) { 371d9f1162dSJed Brown PetscInt mid = lo + (hi - lo)/2; 37259e16d97SJunchao Zhang if (key < X[mid]) hi = mid; 373d9f1162dSJed Brown else lo = mid; 37460e03357SMatthew G Knepley } 37559e16d97SJunchao Zhang *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1); 37660e03357SMatthew G Knepley PetscFunctionReturn(0); 37760e03357SMatthew G Knepley } 37860e03357SMatthew G Knepley 379d2aeb606SJed Brown /*@ 380f1cab4e1SJunchao Zhang PetscCheckDupsInt - Checks if an integer array has duplicates 381f1cab4e1SJunchao Zhang 382f1cab4e1SJunchao Zhang Not Collective 383f1cab4e1SJunchao Zhang 384f1cab4e1SJunchao Zhang Input Parameters: 385f1cab4e1SJunchao Zhang + n - number of values in the array 386f1cab4e1SJunchao Zhang - X - array of integers 387f1cab4e1SJunchao Zhang 388f1cab4e1SJunchao Zhang 389f1cab4e1SJunchao Zhang Output Parameter: 390f1cab4e1SJunchao Zhang . dups - True if the array has dups, otherwise false 391f1cab4e1SJunchao Zhang 392f1cab4e1SJunchao Zhang Level: intermediate 393f1cab4e1SJunchao Zhang 394f1cab4e1SJunchao Zhang .seealso: PetscSortRemoveDupsInt() 395f1cab4e1SJunchao Zhang @*/ 396f1cab4e1SJunchao Zhang PetscErrorCode PetscCheckDupsInt(PetscInt n,const PetscInt X[],PetscBool *dups) 397f1cab4e1SJunchao Zhang { 398f1cab4e1SJunchao Zhang PetscErrorCode ierr; 399f1cab4e1SJunchao Zhang PetscInt i; 400f1cab4e1SJunchao Zhang PetscHSetI ht; 401f1cab4e1SJunchao Zhang PetscBool missing; 402f1cab4e1SJunchao Zhang 403f1cab4e1SJunchao Zhang PetscFunctionBegin; 404f1cab4e1SJunchao Zhang PetscValidPointer(dups,3); 405f1cab4e1SJunchao Zhang *dups = PETSC_FALSE; 406f1cab4e1SJunchao Zhang if (n > 1) { 407f1cab4e1SJunchao Zhang ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 408f1cab4e1SJunchao Zhang ierr = PetscHSetIResize(ht,n);CHKERRQ(ierr); 409f1cab4e1SJunchao Zhang for (i=0; i<n; i++) { 410f1cab4e1SJunchao Zhang ierr = PetscHSetIQueryAdd(ht,X[i],&missing);CHKERRQ(ierr); 411f1cab4e1SJunchao Zhang if (!missing) {*dups = PETSC_TRUE; break;} 412f1cab4e1SJunchao Zhang } 413f1cab4e1SJunchao Zhang ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr); 414f1cab4e1SJunchao Zhang } 415f1cab4e1SJunchao Zhang PetscFunctionReturn(0); 416f1cab4e1SJunchao Zhang } 417f1cab4e1SJunchao Zhang 418f1cab4e1SJunchao Zhang /*@ 419d2aeb606SJed Brown PetscFindMPIInt - Finds MPI integer in a sorted array of integers 420d2aeb606SJed Brown 421d2aeb606SJed Brown Not Collective 422d2aeb606SJed Brown 423d2aeb606SJed Brown Input Parameters: 424d2aeb606SJed Brown + key - the integer to locate 425d2aeb606SJed Brown . n - number of values in the array 42659e16d97SJunchao Zhang - X - array of integers 427d2aeb606SJed Brown 428d2aeb606SJed Brown Output Parameter: 429d2aeb606SJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 430d2aeb606SJed Brown 431d2aeb606SJed Brown Level: intermediate 432d2aeb606SJed Brown 433*676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt() 434d2aeb606SJed Brown @*/ 43559e16d97SJunchao Zhang PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc) 436d2aeb606SJed Brown { 437d2aeb606SJed Brown PetscInt lo = 0,hi = n; 438d2aeb606SJed Brown 439d2aeb606SJed Brown PetscFunctionBegin; 440d2aeb606SJed Brown PetscValidPointer(loc,4); 441d2aeb606SJed Brown if (!n) {*loc = -1; PetscFunctionReturn(0);} 44259e16d97SJunchao Zhang PetscValidPointer(X,3); 4436a7c706bSVaclav Hapla PetscCheckSorted(n,X); 444d2aeb606SJed Brown while (hi - lo > 1) { 445d2aeb606SJed Brown PetscInt mid = lo + (hi - lo)/2; 44659e16d97SJunchao Zhang if (key < X[mid]) hi = mid; 447d2aeb606SJed Brown else lo = mid; 448d2aeb606SJed Brown } 44959e16d97SJunchao Zhang *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1); 450e5c89e4eSSatish Balay PetscFunctionReturn(0); 451e5c89e4eSSatish Balay } 452e5c89e4eSSatish Balay 453e5c89e4eSSatish Balay /*@ 454e5c89e4eSSatish Balay PetscSortIntWithArray - Sorts an array of integers in place in increasing order; 455e5c89e4eSSatish Balay changes a second array to match the sorted first array. 456e5c89e4eSSatish Balay 457e5c89e4eSSatish Balay Not Collective 458e5c89e4eSSatish Balay 459e5c89e4eSSatish Balay Input Parameters: 460e5c89e4eSSatish Balay + n - number of values 46159e16d97SJunchao Zhang . X - array of integers 46259e16d97SJunchao Zhang - Y - second array of integers 463e5c89e4eSSatish Balay 464e5c89e4eSSatish Balay Level: intermediate 465e5c89e4eSSatish Balay 466*676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 467e5c89e4eSSatish Balay @*/ 468b4f531b9SJunchao Zhang PetscErrorCode PetscSortIntWithArray(PetscInt n,PetscInt X[],PetscInt Y[]) 469e5c89e4eSSatish Balay { 470e5c89e4eSSatish Balay PetscErrorCode ierr; 47159e16d97SJunchao Zhang PetscInt pivot,t1,t2; 472e5c89e4eSSatish Balay 473e5c89e4eSSatish Balay PetscFunctionBegin; 47459e16d97SJunchao Zhang QuickSort2(PetscSortIntWithArray,X,Y,n,pivot,t1,t2,ierr); 475c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 476c1f0200aSDmitry Karpeev } 477c1f0200aSDmitry Karpeev 478c1f0200aSDmitry Karpeev /*@ 479c1f0200aSDmitry Karpeev PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order; 480c1f0200aSDmitry Karpeev changes a pair of integer arrays to match the sorted first array. 481c1f0200aSDmitry Karpeev 482c1f0200aSDmitry Karpeev Not Collective 483c1f0200aSDmitry Karpeev 484c1f0200aSDmitry Karpeev Input Parameters: 485c1f0200aSDmitry Karpeev + n - number of values 48659e16d97SJunchao Zhang . X - array of integers 48759e16d97SJunchao Zhang . Y - second array of integers (first array of the pair) 48859e16d97SJunchao Zhang - Z - third array of integers (second array of the pair) 489c1f0200aSDmitry Karpeev 490c1f0200aSDmitry Karpeev Level: intermediate 491c1f0200aSDmitry Karpeev 492*676f2a66SJacob Faibussowitsch .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray(), PetscIntSortSemiOrdered() 493c1f0200aSDmitry Karpeev @*/ 494b4f531b9SJunchao Zhang PetscErrorCode PetscSortIntWithArrayPair(PetscInt n,PetscInt X[],PetscInt Y[],PetscInt Z[]) 495c1f0200aSDmitry Karpeev { 496c1f0200aSDmitry Karpeev PetscErrorCode ierr; 49759e16d97SJunchao Zhang PetscInt pivot,t1,t2,t3; 498c1f0200aSDmitry Karpeev 499c1f0200aSDmitry Karpeev PetscFunctionBegin; 50059e16d97SJunchao Zhang QuickSort3(PetscSortIntWithArrayPair,X,Y,Z,n,pivot,t1,t2,t3,ierr); 501c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 502c1f0200aSDmitry Karpeev } 503c1f0200aSDmitry Karpeev 50417d7d925SStefano Zampini /*@ 5056a7c706bSVaclav Hapla PetscSortedMPIInt - Determines whether the array is sorted. 5066a7c706bSVaclav Hapla 5076a7c706bSVaclav Hapla Not Collective 5086a7c706bSVaclav Hapla 5096a7c706bSVaclav Hapla Input Parameters: 5106a7c706bSVaclav Hapla + n - number of values 5116a7c706bSVaclav Hapla - X - array of integers 5126a7c706bSVaclav Hapla 5136a7c706bSVaclav Hapla Output Parameters: 5146a7c706bSVaclav Hapla . sorted - flag whether the array is sorted 5156a7c706bSVaclav Hapla 5166a7c706bSVaclav Hapla Level: intermediate 5176a7c706bSVaclav Hapla 518*676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortMPIInt(), PetscSortedInt(), PetscSortedReal() 5196a7c706bSVaclav Hapla @*/ 5206a7c706bSVaclav Hapla PetscErrorCode PetscSortedMPIInt(PetscInt n,const PetscMPIInt X[],PetscBool *sorted) 5216a7c706bSVaclav Hapla { 5226a7c706bSVaclav Hapla PetscFunctionBegin; 5236a7c706bSVaclav Hapla PetscSorted(n,X,*sorted); 5246a7c706bSVaclav Hapla PetscFunctionReturn(0); 5256a7c706bSVaclav Hapla } 5266a7c706bSVaclav Hapla 5276a7c706bSVaclav Hapla /*@ 52817d7d925SStefano Zampini PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order. 52917d7d925SStefano Zampini 53017d7d925SStefano Zampini Not Collective 53117d7d925SStefano Zampini 53217d7d925SStefano Zampini Input Parameters: 53317d7d925SStefano Zampini + n - number of values 53459e16d97SJunchao Zhang - X - array of integers 53517d7d925SStefano Zampini 53617d7d925SStefano Zampini Level: intermediate 53717d7d925SStefano Zampini 538*676f2a66SJacob Faibussowitsch Notes: 539*676f2a66SJacob Faibussowitsch This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array 540*676f2a66SJacob Faibussowitsch is completely random. There are exceptions to this and so it is __highly__ recomended that the user benchmark their 541*676f2a66SJacob Faibussowitsch code to see which routine is fastest. 542*676f2a66SJacob Faibussowitsch 543*676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation() 54417d7d925SStefano Zampini @*/ 545b4f531b9SJunchao Zhang PetscErrorCode PetscSortMPIInt(PetscInt n,PetscMPIInt X[]) 54617d7d925SStefano Zampini { 54759e16d97SJunchao Zhang PetscErrorCode ierr; 54859e16d97SJunchao Zhang PetscMPIInt pivot,t1; 54917d7d925SStefano Zampini 55017d7d925SStefano Zampini PetscFunctionBegin; 55159e16d97SJunchao Zhang QuickSort1(PetscSortMPIInt,X,n,pivot,t1,ierr); 55217d7d925SStefano Zampini PetscFunctionReturn(0); 55317d7d925SStefano Zampini } 55417d7d925SStefano Zampini 55517d7d925SStefano Zampini /*@ 55617d7d925SStefano Zampini PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries 55717d7d925SStefano Zampini 55817d7d925SStefano Zampini Not Collective 55917d7d925SStefano Zampini 56017d7d925SStefano Zampini Input Parameters: 56117d7d925SStefano Zampini + n - number of values 56259e16d97SJunchao Zhang - X - array of integers 56317d7d925SStefano Zampini 56417d7d925SStefano Zampini Output Parameter: 56517d7d925SStefano Zampini . n - number of non-redundant values 56617d7d925SStefano Zampini 56717d7d925SStefano Zampini Level: intermediate 56817d7d925SStefano Zampini 56917d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt() 57017d7d925SStefano Zampini @*/ 57159e16d97SJunchao Zhang PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt X[]) 57217d7d925SStefano Zampini { 57317d7d925SStefano Zampini PetscErrorCode ierr; 57417d7d925SStefano Zampini PetscInt i,s = 0,N = *n, b = 0; 57517d7d925SStefano Zampini 57617d7d925SStefano Zampini PetscFunctionBegin; 57759e16d97SJunchao Zhang ierr = PetscSortMPIInt(N,X);CHKERRQ(ierr); 57817d7d925SStefano Zampini for (i=0; i<N-1; i++) { 57959e16d97SJunchao Zhang if (X[b+s+1] != X[b]) { 58059e16d97SJunchao Zhang X[b+1] = X[b+s+1]; b++; 581a297a907SKarl Rupp } else s++; 58217d7d925SStefano Zampini } 58317d7d925SStefano Zampini *n = N - s; 58417d7d925SStefano Zampini PetscFunctionReturn(0); 58517d7d925SStefano Zampini } 586c1f0200aSDmitry Karpeev 5874d615ea0SBarry Smith /*@ 5884d615ea0SBarry Smith PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order; 5894d615ea0SBarry Smith changes a second array to match the sorted first array. 5904d615ea0SBarry Smith 5914d615ea0SBarry Smith Not Collective 5924d615ea0SBarry Smith 5934d615ea0SBarry Smith Input Parameters: 5944d615ea0SBarry Smith + n - number of values 59559e16d97SJunchao Zhang . X - array of integers 59659e16d97SJunchao Zhang - Y - second array of integers 5974d615ea0SBarry Smith 5984d615ea0SBarry Smith Level: intermediate 5994d615ea0SBarry Smith 600*676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 6014d615ea0SBarry Smith @*/ 602b4f531b9SJunchao Zhang PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt X[],PetscMPIInt Y[]) 6034d615ea0SBarry Smith { 6044d615ea0SBarry Smith PetscErrorCode ierr; 60559e16d97SJunchao Zhang PetscMPIInt pivot,t1,t2; 6064d615ea0SBarry Smith 6074d615ea0SBarry Smith PetscFunctionBegin; 60859e16d97SJunchao Zhang QuickSort2(PetscSortMPIIntWithArray,X,Y,n,pivot,t1,t2,ierr); 6095569a785SJunchao Zhang PetscFunctionReturn(0); 6105569a785SJunchao Zhang } 6115569a785SJunchao Zhang 6125569a785SJunchao Zhang /*@ 6135569a785SJunchao Zhang PetscSortMPIIntWithIntArray - Sorts an array of MPI integers in place in increasing order; 6145569a785SJunchao Zhang changes a second array of Petsc intergers to match the sorted first array. 6155569a785SJunchao Zhang 6165569a785SJunchao Zhang Not Collective 6175569a785SJunchao Zhang 6185569a785SJunchao Zhang Input Parameters: 6195569a785SJunchao Zhang + n - number of values 62059e16d97SJunchao Zhang . X - array of MPI integers 62159e16d97SJunchao Zhang - Y - second array of Petsc integers 6225569a785SJunchao Zhang 6235569a785SJunchao Zhang Level: intermediate 6245569a785SJunchao Zhang 6255569a785SJunchao Zhang Notes: this routine is useful when one needs to sort MPI ranks with other integer arrays. 6265569a785SJunchao Zhang 627*676f2a66SJacob Faibussowitsch .seealso: PetscSortMPIIntWithArray(), PetscIntSortSemiOrderedWithArray(), PetscTimSortWithArray() 6285569a785SJunchao Zhang @*/ 629b4f531b9SJunchao Zhang PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n,PetscMPIInt X[],PetscInt Y[]) 6305569a785SJunchao Zhang { 6315569a785SJunchao Zhang PetscErrorCode ierr; 63259e16d97SJunchao Zhang PetscMPIInt pivot,t1; 6335569a785SJunchao Zhang PetscInt t2; 6345569a785SJunchao Zhang 6355569a785SJunchao Zhang PetscFunctionBegin; 63659e16d97SJunchao Zhang QuickSort2(PetscSortMPIIntWithIntArray,X,Y,n,pivot,t1,t2,ierr); 637e5c89e4eSSatish Balay PetscFunctionReturn(0); 638e5c89e4eSSatish Balay } 639e5c89e4eSSatish Balay 640e5c89e4eSSatish Balay /*@ 641e5c89e4eSSatish Balay PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order; 642e5c89e4eSSatish Balay changes a second SCALAR array to match the sorted first INTEGER array. 643e5c89e4eSSatish Balay 644e5c89e4eSSatish Balay Not Collective 645e5c89e4eSSatish Balay 646e5c89e4eSSatish Balay Input Parameters: 647e5c89e4eSSatish Balay + n - number of values 64859e16d97SJunchao Zhang . X - array of integers 64959e16d97SJunchao Zhang - Y - second array of scalars 650e5c89e4eSSatish Balay 651e5c89e4eSSatish Balay Level: intermediate 652e5c89e4eSSatish Balay 653*676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 654e5c89e4eSSatish Balay @*/ 655b4f531b9SJunchao Zhang PetscErrorCode PetscSortIntWithScalarArray(PetscInt n,PetscInt X[],PetscScalar Y[]) 656e5c89e4eSSatish Balay { 657e5c89e4eSSatish Balay PetscErrorCode ierr; 65859e16d97SJunchao Zhang PetscInt pivot,t1; 65959e16d97SJunchao Zhang PetscScalar t2; 660e5c89e4eSSatish Balay 661e5c89e4eSSatish Balay PetscFunctionBegin; 66259e16d97SJunchao Zhang QuickSort2(PetscSortIntWithScalarArray,X,Y,n,pivot,t1,t2,ierr); 66317df18f2SToby Isaac PetscFunctionReturn(0); 66417df18f2SToby Isaac } 66517df18f2SToby Isaac 6666c2863d0SBarry Smith /*@C 66717df18f2SToby Isaac PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order; 66817df18f2SToby Isaac changes a second array to match the sorted first INTEGER array. Unlike other sort routines, the user must 66917df18f2SToby Isaac provide workspace (the size of an element in the data array) to use when sorting. 67017df18f2SToby Isaac 67117df18f2SToby Isaac Not Collective 67217df18f2SToby Isaac 67317df18f2SToby Isaac Input Parameters: 67417df18f2SToby Isaac + n - number of values 67559e16d97SJunchao Zhang . X - array of integers 67659e16d97SJunchao Zhang . Y - second array of data 67717df18f2SToby Isaac . size - sizeof elements in the data array in bytes 67859e16d97SJunchao Zhang - t2 - workspace of "size" bytes used when sorting 67917df18f2SToby Isaac 68017df18f2SToby Isaac Level: intermediate 68117df18f2SToby Isaac 682*676f2a66SJacob Faibussowitsch .seealso: PetscTimSortWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 68317df18f2SToby Isaac @*/ 684b4f531b9SJunchao Zhang PetscErrorCode PetscSortIntWithDataArray(PetscInt n,PetscInt X[],void *Y,size_t size,void *t2) 68517df18f2SToby Isaac { 68617df18f2SToby Isaac PetscErrorCode ierr; 68759e16d97SJunchao Zhang char *YY = (char*)Y; 68859e16d97SJunchao Zhang PetscInt i,j,p,t1,pivot,hi=n-1,l,r; 68917df18f2SToby Isaac 69017df18f2SToby Isaac PetscFunctionBegin; 69117df18f2SToby Isaac if (n<8) { 69259e16d97SJunchao Zhang for (i=0; i<n; i++) { 69359e16d97SJunchao Zhang pivot = X[i]; 69459e16d97SJunchao Zhang for (j=i+1; j<n; j++) { 69559e16d97SJunchao Zhang if (pivot > X[j]) { 69659e16d97SJunchao Zhang SWAP2Data(X[i],X[j],YY+size*i,YY+size*j,t1,t2,size); 69759e16d97SJunchao Zhang pivot = X[i]; 69817df18f2SToby Isaac } 69917df18f2SToby Isaac } 70017df18f2SToby Isaac } 70117df18f2SToby Isaac } else { 70259e16d97SJunchao Zhang /* Two way partition */ 70359e16d97SJunchao Zhang p = MEDIAN(X,hi); 70459e16d97SJunchao Zhang pivot = X[p]; 70559e16d97SJunchao Zhang l = 0; 70659e16d97SJunchao Zhang r = hi; 70759e16d97SJunchao Zhang while (1) { 70859e16d97SJunchao Zhang while (X[l] < pivot) l++; 70959e16d97SJunchao Zhang while (X[r] > pivot) r--; 71059e16d97SJunchao Zhang if (l >= r) {r++; break;} 71159e16d97SJunchao Zhang SWAP2Data(X[l],X[r],YY+size*l,YY+size*r,t1,t2,size); 71259e16d97SJunchao Zhang l++; 71359e16d97SJunchao Zhang r--; 71459e16d97SJunchao Zhang } 71559e16d97SJunchao Zhang ierr = PetscSortIntWithDataArray(l,X,Y,size,t2);CHKERRQ(ierr); 71659e16d97SJunchao Zhang ierr = PetscSortIntWithDataArray(hi-r+1,X+r,YY+size*r,size,t2);CHKERRQ(ierr); 71717df18f2SToby Isaac } 71817df18f2SToby Isaac PetscFunctionReturn(0); 71917df18f2SToby Isaac } 72017df18f2SToby Isaac 72121e72a00SBarry Smith /*@ 72221e72a00SBarry Smith PetscMergeIntArray - Merges two SORTED integer arrays, removes duplicate elements. 72321e72a00SBarry Smith 72421e72a00SBarry Smith Not Collective 72521e72a00SBarry Smith 72621e72a00SBarry Smith Input Parameters: 72721e72a00SBarry Smith + an - number of values in the first array 72821e72a00SBarry Smith . aI - first sorted array of integers 72921e72a00SBarry Smith . bn - number of values in the second array 73021e72a00SBarry Smith - bI - second array of integers 73121e72a00SBarry Smith 73221e72a00SBarry Smith Output Parameters: 73321e72a00SBarry Smith + n - number of values in the merged array 7346c2863d0SBarry Smith - L - merged sorted array, this is allocated if an array is not provided 73521e72a00SBarry Smith 73621e72a00SBarry Smith Level: intermediate 73721e72a00SBarry Smith 73821e72a00SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 73921e72a00SBarry Smith @*/ 7406c2863d0SBarry Smith PetscErrorCode PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L) 74121e72a00SBarry Smith { 74221e72a00SBarry Smith PetscErrorCode ierr; 74321e72a00SBarry Smith PetscInt *L_ = *L, ak, bk, k; 74421e72a00SBarry Smith 74521e72a00SBarry Smith if (!L_) { 74621e72a00SBarry Smith ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr); 74721e72a00SBarry Smith L_ = *L; 74821e72a00SBarry Smith } 74921e72a00SBarry Smith k = ak = bk = 0; 75021e72a00SBarry Smith while (ak < an && bk < bn) { 75121e72a00SBarry Smith if (aI[ak] == bI[bk]) { 75221e72a00SBarry Smith L_[k] = aI[ak]; 75321e72a00SBarry Smith ++ak; 75421e72a00SBarry Smith ++bk; 75521e72a00SBarry Smith ++k; 75621e72a00SBarry Smith } else if (aI[ak] < bI[bk]) { 75721e72a00SBarry Smith L_[k] = aI[ak]; 75821e72a00SBarry Smith ++ak; 75921e72a00SBarry Smith ++k; 76021e72a00SBarry Smith } else { 76121e72a00SBarry Smith L_[k] = bI[bk]; 76221e72a00SBarry Smith ++bk; 76321e72a00SBarry Smith ++k; 76421e72a00SBarry Smith } 76521e72a00SBarry Smith } 76621e72a00SBarry Smith if (ak < an) { 767580bdb30SBarry Smith ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr); 76821e72a00SBarry Smith k += (an-ak); 76921e72a00SBarry Smith } 77021e72a00SBarry Smith if (bk < bn) { 771580bdb30SBarry Smith ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr); 77221e72a00SBarry Smith k += (bn-bk); 77321e72a00SBarry Smith } 77421e72a00SBarry Smith *n = k; 77521e72a00SBarry Smith PetscFunctionReturn(0); 77621e72a00SBarry Smith } 777b4301105SBarry Smith 778c1f0200aSDmitry Karpeev /*@ 77921e72a00SBarry Smith PetscMergeIntArrayPair - Merges two SORTED integer arrays that share NO common values along with an additional array of integers. 780c1f0200aSDmitry Karpeev The additional arrays are the same length as sorted arrays and are merged 781c1f0200aSDmitry Karpeev in the order determined by the merging of the sorted pair. 782c1f0200aSDmitry Karpeev 783c1f0200aSDmitry Karpeev Not Collective 784c1f0200aSDmitry Karpeev 785c1f0200aSDmitry Karpeev Input Parameters: 786c1f0200aSDmitry Karpeev + an - number of values in the first array 787c1f0200aSDmitry Karpeev . aI - first sorted array of integers 788c1f0200aSDmitry Karpeev . aJ - first additional array of integers 789c1f0200aSDmitry Karpeev . bn - number of values in the second array 790c1f0200aSDmitry Karpeev . bI - second array of integers 791c1f0200aSDmitry Karpeev - bJ - second additional of integers 792c1f0200aSDmitry Karpeev 793c1f0200aSDmitry Karpeev Output Parameters: 794c1f0200aSDmitry Karpeev + n - number of values in the merged array (== an + bn) 79514c49006SJed Brown . L - merged sorted array 796c1f0200aSDmitry Karpeev - J - merged additional array 797c1f0200aSDmitry Karpeev 79895452b02SPatrick Sanan Notes: 79995452b02SPatrick Sanan if L or J point to non-null arrays then this routine will assume they are of the approproate size and use them, otherwise this routine will allocate space for them 800c1f0200aSDmitry Karpeev Level: intermediate 801c1f0200aSDmitry Karpeev 802*676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 803c1f0200aSDmitry Karpeev @*/ 8046c2863d0SBarry Smith PetscErrorCode PetscMergeIntArrayPair(PetscInt an,const PetscInt aI[], const PetscInt aJ[], PetscInt bn, const PetscInt bI[], const PetscInt bJ[], PetscInt *n, PetscInt **L, PetscInt **J) 805c1f0200aSDmitry Karpeev { 806c1f0200aSDmitry Karpeev PetscErrorCode ierr; 80770d8d27cSBarry Smith PetscInt n_, *L_, *J_, ak, bk, k; 808c1f0200aSDmitry Karpeev 80914c49006SJed Brown PetscFunctionBegin; 81070d8d27cSBarry Smith PetscValidIntPointer(L,8); 81170d8d27cSBarry Smith PetscValidIntPointer(J,9); 812c1f0200aSDmitry Karpeev n_ = an + bn; 813c1f0200aSDmitry Karpeev *n = n_; 81470d8d27cSBarry Smith if (!*L) { 815785e854fSJed Brown ierr = PetscMalloc1(n_, L);CHKERRQ(ierr); 81670d8d27cSBarry Smith } 817d7aa01a8SHong Zhang L_ = *L; 81870d8d27cSBarry Smith if (!*J) { 81970d8d27cSBarry Smith ierr = PetscMalloc1(n_, J);CHKERRQ(ierr); 820c1f0200aSDmitry Karpeev } 821c1f0200aSDmitry Karpeev J_ = *J; 822c1f0200aSDmitry Karpeev k = ak = bk = 0; 823c1f0200aSDmitry Karpeev while (ak < an && bk < bn) { 824c1f0200aSDmitry Karpeev if (aI[ak] <= bI[bk]) { 825d7aa01a8SHong Zhang L_[k] = aI[ak]; 826c1f0200aSDmitry Karpeev J_[k] = aJ[ak]; 827c1f0200aSDmitry Karpeev ++ak; 828c1f0200aSDmitry Karpeev ++k; 829a297a907SKarl Rupp } else { 830d7aa01a8SHong Zhang L_[k] = bI[bk]; 831c1f0200aSDmitry Karpeev J_[k] = bJ[bk]; 832c1f0200aSDmitry Karpeev ++bk; 833c1f0200aSDmitry Karpeev ++k; 834c1f0200aSDmitry Karpeev } 835c1f0200aSDmitry Karpeev } 836c1f0200aSDmitry Karpeev if (ak < an) { 837580bdb30SBarry Smith ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr); 838580bdb30SBarry Smith ierr = PetscArraycpy(J_+k,aJ+ak,an-ak);CHKERRQ(ierr); 839c1f0200aSDmitry Karpeev k += (an-ak); 840c1f0200aSDmitry Karpeev } 841c1f0200aSDmitry Karpeev if (bk < bn) { 842580bdb30SBarry Smith ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr); 843580bdb30SBarry Smith ierr = PetscArraycpy(J_+k,bJ+bk,bn-bk);CHKERRQ(ierr); 844c1f0200aSDmitry Karpeev } 845c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 846c1f0200aSDmitry Karpeev } 847c1f0200aSDmitry Karpeev 848e498c390SJed Brown /*@ 849e498c390SJed Brown PetscMergeMPIIntArray - Merges two SORTED integer arrays. 850e498c390SJed Brown 851e498c390SJed Brown Not Collective 852e498c390SJed Brown 853e498c390SJed Brown Input Parameters: 854e498c390SJed Brown + an - number of values in the first array 855e498c390SJed Brown . aI - first sorted array of integers 856e498c390SJed Brown . bn - number of values in the second array 857e498c390SJed Brown - bI - second array of integers 858e498c390SJed Brown 859e498c390SJed Brown Output Parameters: 860e498c390SJed Brown + n - number of values in the merged array (<= an + bn) 861e498c390SJed Brown - L - merged sorted array, allocated if address of NULL pointer is passed 862e498c390SJed Brown 863e498c390SJed Brown Level: intermediate 864e498c390SJed Brown 865*676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 866e498c390SJed Brown @*/ 867e498c390SJed Brown PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L) 868e498c390SJed Brown { 869e498c390SJed Brown PetscErrorCode ierr; 870e498c390SJed Brown PetscInt ai,bi,k; 871e498c390SJed Brown 872e498c390SJed Brown PetscFunctionBegin; 873e498c390SJed Brown if (!*L) {ierr = PetscMalloc1((an+bn),L);CHKERRQ(ierr);} 874e498c390SJed Brown for (ai=0,bi=0,k=0; ai<an || bi<bn;) { 875e498c390SJed Brown PetscInt t = -1; 876e498c390SJed Brown for (; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai]; 877e498c390SJed Brown for (; bi<bn && bI[bi] == t; bi++); 878e498c390SJed Brown for (; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi]; 879e498c390SJed Brown for (; ai<an && aI[ai] == t; ai++); 880e498c390SJed Brown } 881e498c390SJed Brown *n = k; 882e498c390SJed Brown PetscFunctionReturn(0); 883e498c390SJed Brown } 884e5c89e4eSSatish Balay 8856c2863d0SBarry Smith /*@C 88642eaa7f3SBarry Smith PetscProcessTree - Prepares tree data to be displayed graphically 88742eaa7f3SBarry Smith 88842eaa7f3SBarry Smith Not Collective 88942eaa7f3SBarry Smith 89042eaa7f3SBarry Smith Input Parameters: 89142eaa7f3SBarry Smith + n - number of values 89242eaa7f3SBarry Smith . mask - indicates those entries in the tree, location 0 is always masked 89342eaa7f3SBarry Smith - parentid - indicates the parent of each entry 89442eaa7f3SBarry Smith 89542eaa7f3SBarry Smith Output Parameters: 89642eaa7f3SBarry Smith + Nlevels - the number of levels 89742eaa7f3SBarry Smith . Level - for each node tells its level 89842eaa7f3SBarry Smith . Levelcnts - the number of nodes on each level 89942eaa7f3SBarry Smith . Idbylevel - a list of ids on each of the levels, first level followed by second etc 90042eaa7f3SBarry Smith - Column - for each id tells its column index 90142eaa7f3SBarry Smith 9026c2863d0SBarry Smith Level: developer 90342eaa7f3SBarry Smith 90495452b02SPatrick Sanan Notes: 90595452b02SPatrick Sanan This code is not currently used 90642eaa7f3SBarry Smith 90742eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation() 90842eaa7f3SBarry Smith @*/ 9097087cfbeSBarry Smith PetscErrorCode PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column) 91042eaa7f3SBarry Smith { 91142eaa7f3SBarry Smith PetscInt i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column; 91242eaa7f3SBarry Smith PetscErrorCode ierr; 913ace3abfcSBarry Smith PetscBool done = PETSC_FALSE; 91442eaa7f3SBarry Smith 91542eaa7f3SBarry Smith PetscFunctionBegin; 91642eaa7f3SBarry Smith if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set"); 91742eaa7f3SBarry Smith for (i=0; i<n; i++) { 91842eaa7f3SBarry Smith if (mask[i]) continue; 91942eaa7f3SBarry Smith if (parentid[i] == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent"); 92042eaa7f3SBarry Smith if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked"); 92142eaa7f3SBarry Smith } 92242eaa7f3SBarry Smith 92342eaa7f3SBarry Smith for (i=0; i<n; i++) { 92442eaa7f3SBarry Smith if (!mask[i]) nmask++; 92542eaa7f3SBarry Smith } 92642eaa7f3SBarry Smith 92742eaa7f3SBarry Smith /* determine the level in the tree of each node */ 9281795a4d1SJed Brown ierr = PetscCalloc1(n,&level);CHKERRQ(ierr); 929a297a907SKarl Rupp 93042eaa7f3SBarry Smith level[0] = 1; 93142eaa7f3SBarry Smith while (!done) { 93242eaa7f3SBarry Smith done = PETSC_TRUE; 93342eaa7f3SBarry Smith for (i=0; i<n; i++) { 93442eaa7f3SBarry Smith if (mask[i]) continue; 93542eaa7f3SBarry Smith if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1; 93642eaa7f3SBarry Smith else if (!level[i]) done = PETSC_FALSE; 93742eaa7f3SBarry Smith } 93842eaa7f3SBarry Smith } 93942eaa7f3SBarry Smith for (i=0; i<n; i++) { 94042eaa7f3SBarry Smith level[i]--; 94142eaa7f3SBarry Smith nlevels = PetscMax(nlevels,level[i]); 94242eaa7f3SBarry Smith } 94342eaa7f3SBarry Smith 94442eaa7f3SBarry Smith /* count the number of nodes on each level and its max */ 9451795a4d1SJed Brown ierr = PetscCalloc1(nlevels,&levelcnt);CHKERRQ(ierr); 94642eaa7f3SBarry Smith for (i=0; i<n; i++) { 94742eaa7f3SBarry Smith if (mask[i]) continue; 94842eaa7f3SBarry Smith levelcnt[level[i]-1]++; 94942eaa7f3SBarry Smith } 950a297a907SKarl Rupp for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]); 95142eaa7f3SBarry Smith 95242eaa7f3SBarry Smith /* for each level sort the ids by the parent id */ 953dcca6d9dSJed Brown ierr = PetscMalloc2(levelmax,&workid,levelmax,&workparentid);CHKERRQ(ierr); 954785e854fSJed Brown ierr = PetscMalloc1(nmask,&idbylevel);CHKERRQ(ierr); 95542eaa7f3SBarry Smith for (j=1; j<=nlevels;j++) { 95642eaa7f3SBarry Smith cnt = 0; 95742eaa7f3SBarry Smith for (i=0; i<n; i++) { 95842eaa7f3SBarry Smith if (mask[i]) continue; 95942eaa7f3SBarry Smith if (level[i] != j) continue; 96042eaa7f3SBarry Smith workid[cnt] = i; 96142eaa7f3SBarry Smith workparentid[cnt++] = parentid[i]; 96242eaa7f3SBarry Smith } 96342eaa7f3SBarry Smith /* PetscIntView(cnt,workparentid,0); 96442eaa7f3SBarry Smith PetscIntView(cnt,workid,0); 96542eaa7f3SBarry Smith ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr); 96642eaa7f3SBarry Smith PetscIntView(cnt,workparentid,0); 96742eaa7f3SBarry Smith PetscIntView(cnt,workid,0);*/ 968580bdb30SBarry Smith ierr = PetscArraycpy(idbylevel+tcnt,workid,cnt);CHKERRQ(ierr); 96942eaa7f3SBarry Smith tcnt += cnt; 97042eaa7f3SBarry Smith } 97142eaa7f3SBarry Smith if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes"); 97242eaa7f3SBarry Smith ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr); 97342eaa7f3SBarry Smith 97442eaa7f3SBarry Smith /* for each node list its column */ 975785e854fSJed Brown ierr = PetscMalloc1(n,&column);CHKERRQ(ierr); 97642eaa7f3SBarry Smith cnt = 0; 97742eaa7f3SBarry Smith for (j=0; j<nlevels; j++) { 97842eaa7f3SBarry Smith for (i=0; i<levelcnt[j]; i++) { 97942eaa7f3SBarry Smith column[idbylevel[cnt++]] = i; 98042eaa7f3SBarry Smith } 98142eaa7f3SBarry Smith } 98242eaa7f3SBarry Smith 98342eaa7f3SBarry Smith *Nlevels = nlevels; 98442eaa7f3SBarry Smith *Level = level; 98542eaa7f3SBarry Smith *Levelcnt = levelcnt; 98642eaa7f3SBarry Smith *Idbylevel = idbylevel; 98742eaa7f3SBarry Smith *Column = column; 98842eaa7f3SBarry Smith PetscFunctionReturn(0); 98942eaa7f3SBarry Smith } 990ce605777SToby Isaac 991ce605777SToby Isaac /*@ 992ce605777SToby Isaac PetscParallelSortedInt - Check whether an integer array, distributed over a communicator, is globally sorted. 993ce605777SToby Isaac 994ce605777SToby Isaac Collective 995ce605777SToby Isaac 996ce605777SToby Isaac Input Parameters: 997ce605777SToby Isaac + comm - the MPI communicator 998ce605777SToby Isaac . n - the local number of integers 999ce605777SToby Isaac - keys - the local array of integers 1000ce605777SToby Isaac 1001ce605777SToby Isaac Output Parameters: 1002ce605777SToby Isaac . is_sorted - whether the array is globally sorted 1003ce605777SToby Isaac 1004ce605777SToby Isaac Level: developer 1005ce605777SToby Isaac 1006ce605777SToby Isaac .seealso: PetscParallelSortInt() 1007ce605777SToby Isaac @*/ 1008ce605777SToby Isaac PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted) 1009ce605777SToby Isaac { 1010ce605777SToby Isaac PetscBool sorted; 1011ce605777SToby Isaac PetscInt i, min, max, prevmax; 10122a1da528SToby Isaac PetscMPIInt rank; 1013ce605777SToby Isaac PetscErrorCode ierr; 1014ce605777SToby Isaac 1015ce605777SToby Isaac PetscFunctionBegin; 1016ce605777SToby Isaac sorted = PETSC_TRUE; 1017ce605777SToby Isaac min = PETSC_MAX_INT; 1018ce605777SToby Isaac max = PETSC_MIN_INT; 1019ce605777SToby Isaac if (n) { 1020ce605777SToby Isaac min = keys[0]; 1021ce605777SToby Isaac max = keys[0]; 1022ce605777SToby Isaac } 1023ce605777SToby Isaac for (i = 1; i < n; i++) { 1024ce605777SToby Isaac if (keys[i] < keys[i - 1]) break; 1025ce605777SToby Isaac min = PetscMin(min,keys[i]); 1026ce605777SToby Isaac max = PetscMax(max,keys[i]); 1027ce605777SToby Isaac } 1028ce605777SToby Isaac if (i < n) sorted = PETSC_FALSE; 1029ce605777SToby Isaac prevmax = PETSC_MIN_INT; 1030ce605777SToby Isaac ierr = MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 10312a1da528SToby Isaac ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 10322a1da528SToby Isaac if (!rank) prevmax = PETSC_MIN_INT; 1033ce605777SToby Isaac if (prevmax > min) sorted = PETSC_FALSE; 1034ce605777SToby Isaac ierr = MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr); 1035ce605777SToby Isaac PetscFunctionReturn(0); 1036ce605777SToby Isaac } 1037