/* This file contains routines for sorting integers. Values are sorted in place. One can use src/sys/tests/ex52.c for benchmarking. */ #include /*I "petscsys.h" I*/ #include #define MEDIAN3(v,a,b,c) \ (v[a] */ #define SWAP2Data(a,b,c,d,t1,t2,siz) \ do { \ PetscErrorCode ierr; \ t1=a; a=b; b=t1; \ ierr = PetscMemcpy(t2,c,siz);CHKERRQ(ierr); \ ierr = PetscMemcpy(c,d,siz);CHKERRQ(ierr); \ ierr = PetscMemcpy(d,t2,siz);CHKERRQ(ierr); \ } while (0) /* Partition X[lo,hi] into two parts: X[lo,l) <= pivot; X[r,hi] > pivot Input Parameters: + X - array to partition . pivot - a pivot of X[] . t1 - temp variable for X - lo,hi - lower and upper bound of the array Output Parameters: + l,r - of type PetscInt Notes: The TwoWayPartition2/3 variants also partition other arrays along with X. These arrays can have different types, so they provide their own temp t2,t3 */ #define TwoWayPartition1(X,pivot,t1,lo,hi,l,r) \ do { \ l = lo; \ r = hi; \ while (1) { \ while (X[l] < pivot) l++; \ while (X[r] > pivot) r--; \ if (l >= r) {r++; break;} \ SWAP1(X[l],X[r],t1); \ l++; \ r--; \ } \ } while (0) /* Partition X[lo,hi] into two parts: X[lo,l) >= pivot; X[r,hi] < pivot Input Parameters: + X - array to partition . pivot - a pivot of X[] . t1 - temp variable for X - lo,hi - lower and upper bound of the array Output Parameters: + l,r - of type PetscInt Notes: The TwoWayPartition2/3 variants also partition other arrays along with X. These arrays can have different types, so they provide their own temp t2,t3 */ #define TwoWayPartitionReverse1(X,pivot,t1,lo,hi,l,r) \ do { \ l = lo; \ r = hi; \ while (1) { \ while (X[l] > pivot) l++; \ while (X[r] < pivot) r--; \ if (l >= r) {r++; break;} \ SWAP1(X[l],X[r],t1); \ l++; \ r--; \ } \ } while (0) #define TwoWayPartition2(X,Y,pivot,t1,t2,lo,hi,l,r) \ do { \ l = lo; \ r = hi; \ while (1) { \ while (X[l] < pivot) l++; \ while (X[r] > pivot) r--; \ if (l >= r) {r++; break;} \ SWAP2(X[l],X[r],Y[l],Y[r],t1,t2); \ l++; \ r--; \ } \ } while (0) #define TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,lo,hi,l,r) \ do { \ l = lo; \ r = hi; \ while (1) { \ while (X[l] < pivot) l++; \ while (X[r] > pivot) r--; \ if (l >= r) {r++; break;} \ SWAP3(X[l],X[r],Y[l],Y[r],Z[l],Z[r],t1,t2,t3); \ l++; \ r--; \ } \ } while (0) /* Templates for similar functions used below */ #define QuickSort1(FuncName,X,n,pivot,t1,ierr) \ do { \ PetscInt i,j,p,l,r,hi=n-1; \ if (n < 8) { \ for (i=0; i X[j]) { \ SWAP1(X[i],X[j],t1); \ pivot = X[i]; \ } \ } \ } \ } else { \ p = MEDIAN(X,hi); \ pivot = X[p]; \ TwoWayPartition1(X,pivot,t1,0,hi,l,r); \ ierr = FuncName(l,X);CHKERRQ(ierr); \ ierr = FuncName(hi-r+1,X+r);CHKERRQ(ierr); \ } \ } while (0) /* Templates for similar functions used below */ #define QuickSortReverse1(FuncName,X,n,pivot,t1,ierr) \ do { \ PetscInt i,j,p,l,r,hi=n-1; \ if (n < 8) { \ for (i=0; i X[j]) { \ SWAP2(X[i],X[j],Y[i],Y[j],t1,t2); \ pivot = X[i]; \ } \ } \ } \ } else { \ p = MEDIAN(X,hi); \ pivot = X[p]; \ TwoWayPartition2(X,Y,pivot,t1,t2,0,hi,l,r); \ ierr = FuncName(l,X,Y);CHKERRQ(ierr); \ ierr = FuncName(hi-r+1,X+r,Y+r);CHKERRQ(ierr); \ } \ } while (0) #define QuickSort3(FuncName,X,Y,Z,n,pivot,t1,t2,t3,ierr) \ do { \ PetscInt i,j,p,l,r,hi=n-1; \ if (n < 8) { \ for (i=0; i X[j]) { \ SWAP3(X[i],X[j],Y[i],Y[j],Z[i],Z[j],t1,t2,t3); \ pivot = X[i]; \ } \ } \ } \ } else { \ p = MEDIAN(X,hi); \ pivot = X[p]; \ TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,0,hi,l,r); \ ierr = FuncName(l,X,Y,Z);CHKERRQ(ierr); \ ierr = FuncName(hi-r+1,X+r,Y+r,Z+r);CHKERRQ(ierr); \ } \ } while (0) /*@ PetscSortedInt - Determines whether the array is sorted. Not Collective Input Parameters: + n - number of values - X - array of integers Output Parameters: . sorted - flag whether the array is sorted Level: intermediate .seealso: PetscSortInt(), PetscSortedMPIInt(), PetscSortedReal() @*/ PetscErrorCode PetscSortedInt(PetscInt n,const PetscInt X[],PetscBool *sorted) { PetscFunctionBegin; PetscSorted(n,X,*sorted); PetscFunctionReturn(0); } /*@ PetscSortInt - Sorts an array of integers in place in increasing order. Not Collective Input Parameters: + n - number of values - X - array of integers Notes: This function serves as an alternative to PetscIntSortSemiOrdered(), and may perform faster especially if the array is completely random. There are exceptions to this and so it is __highly__ recomended that the user benchmark their code to see which routine is fastest. Level: intermediate .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation() @*/ PetscErrorCode PetscSortInt(PetscInt n,PetscInt X[]) { PetscErrorCode ierr; PetscInt pivot,t1; PetscFunctionBegin; QuickSort1(PetscSortInt,X,n,pivot,t1,ierr); PetscFunctionReturn(0); } /*@ PetscSortReverseInt - Sorts an array of integers in place in decreasing order. Not Collective Input Parameters: + n - number of values - X - array of integers Level: intermediate .seealso: PetscIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithPermutation() @*/ PetscErrorCode PetscSortReverseInt(PetscInt n,PetscInt X[]) { PetscErrorCode ierr; PetscInt pivot,t1; PetscFunctionBegin; QuickSortReverse1(PetscSortReverseInt,X,n,pivot,t1,ierr); PetscFunctionReturn(0); } /*@ PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted input array Not Collective Input Parameters: + n - number of values - X - sorted array of integers Output Parameter: . n - number of non-redundant values Level: intermediate .seealso: PetscSortInt() @*/ PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n,PetscInt X[]) { PetscInt i,s = 0,N = *n, b = 0; PetscFunctionBegin; PetscCheckSorted(*n,X); for (i=0; i 1) { PetscInt mid = lo + (hi - lo)/2; if (key < X[mid]) hi = mid; else lo = mid; } *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1); PetscFunctionReturn(0); } /*@ PetscCheckDupsInt - Checks if an integer array has duplicates Not Collective Input Parameters: + n - number of values in the array - X - array of integers Output Parameter: . dups - True if the array has dups, otherwise false Level: intermediate .seealso: PetscSortRemoveDupsInt() @*/ PetscErrorCode PetscCheckDupsInt(PetscInt n,const PetscInt X[],PetscBool *dups) { PetscErrorCode ierr; PetscInt i; PetscHSetI ht; PetscBool missing; PetscFunctionBegin; PetscValidPointer(dups,3); *dups = PETSC_FALSE; if (n > 1) { ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); ierr = PetscHSetIResize(ht,n);CHKERRQ(ierr); for (i=0; i 1) { PetscInt mid = lo + (hi - lo)/2; if (key < X[mid]) hi = mid; else lo = mid; } *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1); PetscFunctionReturn(0); } /*@ PetscSortIntWithArray - Sorts an array of integers in place in increasing order; changes a second array to match the sorted first array. Not Collective Input Parameters: + n - number of values . X - array of integers - Y - second array of integers Level: intermediate .seealso: PetscIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() @*/ PetscErrorCode PetscSortIntWithArray(PetscInt n,PetscInt X[],PetscInt Y[]) { PetscErrorCode ierr; PetscInt pivot,t1,t2; PetscFunctionBegin; QuickSort2(PetscSortIntWithArray,X,Y,n,pivot,t1,t2,ierr); PetscFunctionReturn(0); } /*@ PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order; changes a pair of integer arrays to match the sorted first array. Not Collective Input Parameters: + n - number of values . X - array of integers . Y - second array of integers (first array of the pair) - Z - third array of integers (second array of the pair) Level: intermediate .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray(), PetscIntSortSemiOrdered() @*/ PetscErrorCode PetscSortIntWithArrayPair(PetscInt n,PetscInt X[],PetscInt Y[],PetscInt Z[]) { PetscErrorCode ierr; PetscInt pivot,t1,t2,t3; PetscFunctionBegin; QuickSort3(PetscSortIntWithArrayPair,X,Y,Z,n,pivot,t1,t2,t3,ierr); PetscFunctionReturn(0); } /*@ PetscSortedMPIInt - Determines whether the array is sorted. Not Collective Input Parameters: + n - number of values - X - array of integers Output Parameters: . sorted - flag whether the array is sorted Level: intermediate .seealso: PetscMPIIntSortSemiOrdered(), PetscSortMPIInt(), PetscSortedInt(), PetscSortedReal() @*/ PetscErrorCode PetscSortedMPIInt(PetscInt n,const PetscMPIInt X[],PetscBool *sorted) { PetscFunctionBegin; PetscSorted(n,X,*sorted); PetscFunctionReturn(0); } /*@ PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order. Not Collective Input Parameters: + n - number of values - X - array of integers Level: intermediate Notes: This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array is completely random. There are exceptions to this and so it is __highly__ recomended that the user benchmark their code to see which routine is fastest. .seealso: PetscMPIIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation() @*/ PetscErrorCode PetscSortMPIInt(PetscInt n,PetscMPIInt X[]) { PetscErrorCode ierr; PetscMPIInt pivot,t1; PetscFunctionBegin; QuickSort1(PetscSortMPIInt,X,n,pivot,t1,ierr); PetscFunctionReturn(0); } /*@ PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries Not Collective Input Parameters: + n - number of values - X - array of integers Output Parameter: . n - number of non-redundant values Level: intermediate .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt() @*/ PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt X[]) { PetscErrorCode ierr; PetscInt i,s = 0,N = *n, b = 0; PetscFunctionBegin; ierr = PetscSortMPIInt(N,X);CHKERRQ(ierr); for (i=0; i X[j]) { SWAP2Data(X[i],X[j],YY+size*i,YY+size*j,t1,t2,size); pivot = X[i]; } } } } else { /* Two way partition */ p = MEDIAN(X,hi); pivot = X[p]; l = 0; r = hi; while (1) { while (X[l] < pivot) l++; while (X[r] > pivot) r--; if (l >= r) {r++; break;} SWAP2Data(X[l],X[r],YY+size*l,YY+size*r,t1,t2,size); l++; r--; } ierr = PetscSortIntWithDataArray(l,X,Y,size,t2);CHKERRQ(ierr); ierr = PetscSortIntWithDataArray(hi-r+1,X+r,YY+size*r,size,t2);CHKERRQ(ierr); } PetscFunctionReturn(0); } /*@ PetscMergeIntArray - Merges two SORTED integer arrays, removes duplicate elements. Not Collective Input Parameters: + an - number of values in the first array . aI - first sorted array of integers . bn - number of values in the second array - bI - second array of integers Output Parameters: + n - number of values in the merged array - L - merged sorted array, this is allocated if an array is not provided Level: intermediate .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() @*/ PetscErrorCode PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L) { PetscErrorCode ierr; PetscInt *L_ = *L, ak, bk, k; if (!L_) { ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr); L_ = *L; } k = ak = bk = 0; while (ak < an && bk < bn) { if (aI[ak] == bI[bk]) { L_[k] = aI[ak]; ++ak; ++bk; ++k; } else if (aI[ak] < bI[bk]) { L_[k] = aI[ak]; ++ak; ++k; } else { L_[k] = bI[bk]; ++bk; ++k; } } if (ak < an) { ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr); k += (an-ak); } if (bk < bn) { ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr); k += (bn-bk); } *n = k; PetscFunctionReturn(0); } /*@ PetscMergeIntArrayPair - Merges two SORTED integer arrays that share NO common values along with an additional array of integers. The additional arrays are the same length as sorted arrays and are merged in the order determined by the merging of the sorted pair. Not Collective Input Parameters: + an - number of values in the first array . aI - first sorted array of integers . aJ - first additional array of integers . bn - number of values in the second array . bI - second array of integers - bJ - second additional of integers Output Parameters: + n - number of values in the merged array (== an + bn) . L - merged sorted array - J - merged additional array Notes: 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 Level: intermediate .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() @*/ PetscErrorCode PetscMergeIntArrayPair(PetscInt an,const PetscInt aI[], const PetscInt aJ[], PetscInt bn, const PetscInt bI[], const PetscInt bJ[], PetscInt *n, PetscInt **L, PetscInt **J) { PetscErrorCode ierr; PetscInt n_, *L_, *J_, ak, bk, k; PetscFunctionBegin; PetscValidIntPointer(L,8); PetscValidIntPointer(J,9); n_ = an + bn; *n = n_; if (!*L) { ierr = PetscMalloc1(n_, L);CHKERRQ(ierr); } L_ = *L; if (!*J) { ierr = PetscMalloc1(n_, J);CHKERRQ(ierr); } J_ = *J; k = ak = bk = 0; while (ak < an && bk < bn) { if (aI[ak] <= bI[bk]) { L_[k] = aI[ak]; J_[k] = aJ[ak]; ++ak; ++k; } else { L_[k] = bI[bk]; J_[k] = bJ[bk]; ++bk; ++k; } } if (ak < an) { ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr); ierr = PetscArraycpy(J_+k,aJ+ak,an-ak);CHKERRQ(ierr); k += (an-ak); } if (bk < bn) { ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr); ierr = PetscArraycpy(J_+k,bJ+bk,bn-bk);CHKERRQ(ierr); } PetscFunctionReturn(0); } /*@ PetscMergeMPIIntArray - Merges two SORTED integer arrays. Not Collective Input Parameters: + an - number of values in the first array . aI - first sorted array of integers . bn - number of values in the second array - bI - second array of integers Output Parameters: + n - number of values in the merged array (<= an + bn) - L - merged sorted array, allocated if address of NULL pointer is passed Level: intermediate .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() @*/ PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L) { PetscErrorCode ierr; PetscInt ai,bi,k; PetscFunctionBegin; if (!*L) {ierr = PetscMalloc1((an+bn),L);CHKERRQ(ierr);} for (ai=0,bi=0,k=0; ai min) sorted = PETSC_FALSE; ierr = MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr); PetscFunctionReturn(0); }