17d0a6c19SBarry Smith 2e5c89e4eSSatish Balay /* 3e5c89e4eSSatish Balay This file contains routines for sorting integers. Values are sorted in place. 4e5c89e4eSSatish Balay */ 5af0996ceSBarry Smith #include <petsc/private/petscimpl.h> /*I "petscsys.h" I*/ 6e5c89e4eSSatish Balay 7e5c89e4eSSatish Balay #define SWAP(a,b,t) {t=a;a=b;b=t;} 8e5c89e4eSSatish Balay 9a8582168SJed Brown #define MEDIAN3(v,a,b,c) \ 10a8582168SJed Brown (v[a]<v[b] \ 11a8582168SJed Brown ? (v[b]<v[c] \ 12a8582168SJed Brown ? b \ 13a8582168SJed Brown : (v[a]<v[c] ? c : a)) \ 14a8582168SJed Brown : (v[c]<v[b] \ 15a8582168SJed Brown ? b \ 16a8582168SJed Brown : (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 20e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 21e5c89e4eSSatish Balay 22e5c89e4eSSatish Balay /* 23e5c89e4eSSatish Balay A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 24e5c89e4eSSatish Balay Assumes 0 origin for v, number of elements = right+1 (right is index of 25e5c89e4eSSatish Balay right-most member). 26e5c89e4eSSatish Balay */ 27ef8e3583SJed Brown static void PetscSortInt_Private(PetscInt *v,PetscInt right) 28e5c89e4eSSatish Balay { 29ef8e3583SJed Brown PetscInt i,j,pivot,tmp; 30e5c89e4eSSatish Balay 31e5c89e4eSSatish Balay if (right <= 1) { 32e5c89e4eSSatish Balay if (right == 1) { 33e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP(v[0],v[1],tmp); 34e5c89e4eSSatish Balay } 35ef8e3583SJed Brown return; 36e5c89e4eSSatish Balay } 37ef8e3583SJed Brown i = MEDIAN(v,right); /* Choose a pivot */ 38ef8e3583SJed Brown SWAP(v[0],v[i],tmp); /* Move it out of the way */ 39ef8e3583SJed Brown pivot = v[0]; 40ef8e3583SJed Brown for (i=0,j=right+1;; ) { 41ef8e3583SJed Brown while (++i < j && v[i] <= pivot) ; /* Scan from the left */ 42ef8e3583SJed Brown while (v[--j] > pivot) ; /* Scan from the right */ 43ef8e3583SJed Brown if (i >= j) break; 44ef8e3583SJed Brown SWAP(v[i],v[j],tmp); 45e5c89e4eSSatish Balay } 46ef8e3583SJed Brown SWAP(v[0],v[j],tmp); /* Put pivot back in place. */ 47ef8e3583SJed Brown PetscSortInt_Private(v,j-1); 48ef8e3583SJed Brown PetscSortInt_Private(v+j+1,right-(j+1)); 49e5c89e4eSSatish Balay } 50e5c89e4eSSatish Balay 51e5c89e4eSSatish Balay /*@ 52e5c89e4eSSatish Balay PetscSortInt - Sorts an array of integers in place in increasing order. 53e5c89e4eSSatish Balay 54e5c89e4eSSatish Balay Not Collective 55e5c89e4eSSatish Balay 56e5c89e4eSSatish Balay Input Parameters: 57e5c89e4eSSatish Balay + n - number of values 58e5c89e4eSSatish Balay - i - array of integers 59e5c89e4eSSatish Balay 60e5c89e4eSSatish Balay Level: intermediate 61e5c89e4eSSatish Balay 62e5c89e4eSSatish Balay Concepts: sorting^ints 63e5c89e4eSSatish Balay 64e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntWithPermutation() 65e5c89e4eSSatish Balay @*/ 667087cfbeSBarry Smith PetscErrorCode PetscSortInt(PetscInt n,PetscInt i[]) 67e5c89e4eSSatish Balay { 68e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 69e5c89e4eSSatish Balay 70e5c89e4eSSatish Balay PetscFunctionBegin; 71e5c89e4eSSatish Balay if (n<8) { 72e5c89e4eSSatish Balay for (k=0; k<n; k++) { 73e5c89e4eSSatish Balay ik = i[k]; 74e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 75e5c89e4eSSatish Balay if (ik > i[j]) { 76e5c89e4eSSatish Balay SWAP(i[k],i[j],tmp); 77e5c89e4eSSatish Balay ik = i[k]; 78e5c89e4eSSatish Balay } 79e5c89e4eSSatish Balay } 80e5c89e4eSSatish Balay } 81a297a907SKarl Rupp } else PetscSortInt_Private(i,n-1); 82e5c89e4eSSatish Balay PetscFunctionReturn(0); 83e5c89e4eSSatish Balay } 84e5c89e4eSSatish Balay 851db36a52SBarry Smith /*@ 8622ab5688SLisandro Dalcin PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted input array 8722ab5688SLisandro Dalcin 8822ab5688SLisandro Dalcin Not Collective 8922ab5688SLisandro Dalcin 9022ab5688SLisandro Dalcin Input Parameters: 9122ab5688SLisandro Dalcin + n - number of values 9222ab5688SLisandro Dalcin - ii - sorted array of integers 9322ab5688SLisandro Dalcin 9422ab5688SLisandro Dalcin Output Parameter: 9522ab5688SLisandro Dalcin . n - number of non-redundant values 9622ab5688SLisandro Dalcin 9722ab5688SLisandro Dalcin Level: intermediate 9822ab5688SLisandro Dalcin 9922ab5688SLisandro Dalcin Concepts: sorting^ints 10022ab5688SLisandro Dalcin 10122ab5688SLisandro Dalcin .seealso: PetscSortInt() 10222ab5688SLisandro Dalcin @*/ 10322ab5688SLisandro Dalcin PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n,PetscInt ii[]) 10422ab5688SLisandro Dalcin { 10522ab5688SLisandro Dalcin PetscInt i,s = 0,N = *n, b = 0; 10622ab5688SLisandro Dalcin 10722ab5688SLisandro Dalcin PetscFunctionBegin; 10822ab5688SLisandro Dalcin for (i=0; i<N-1; i++) { 10922ab5688SLisandro Dalcin if (PetscUnlikely(ii[b+s+1] < ii[b])) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Input array is not sorted"); 11022ab5688SLisandro Dalcin if (ii[b+s+1] != ii[b]) { 11122ab5688SLisandro Dalcin ii[b+1] = ii[b+s+1]; b++; 11222ab5688SLisandro Dalcin } else s++; 11322ab5688SLisandro Dalcin } 11422ab5688SLisandro Dalcin *n = N - s; 11522ab5688SLisandro Dalcin PetscFunctionReturn(0); 11622ab5688SLisandro Dalcin } 11722ab5688SLisandro Dalcin 11822ab5688SLisandro Dalcin /*@ 1191db36a52SBarry Smith PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries 1201db36a52SBarry Smith 1211db36a52SBarry Smith Not Collective 1221db36a52SBarry Smith 1231db36a52SBarry Smith Input Parameters: 1241db36a52SBarry Smith + n - number of values 1251db36a52SBarry Smith - ii - array of integers 1261db36a52SBarry Smith 1271db36a52SBarry Smith Output Parameter: 1281db36a52SBarry Smith . n - number of non-redundant values 1291db36a52SBarry Smith 1301db36a52SBarry Smith Level: intermediate 1311db36a52SBarry Smith 1321db36a52SBarry Smith Concepts: sorting^ints 1331db36a52SBarry Smith 13422ab5688SLisandro Dalcin .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortedRemoveDupsInt() 1351db36a52SBarry Smith @*/ 1367087cfbeSBarry Smith PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n,PetscInt ii[]) 1371db36a52SBarry Smith { 1381db36a52SBarry Smith PetscErrorCode ierr; 1391db36a52SBarry Smith 1401db36a52SBarry Smith PetscFunctionBegin; 14122ab5688SLisandro Dalcin ierr = PetscSortInt(*n,ii);CHKERRQ(ierr); 14222ab5688SLisandro Dalcin ierr = PetscSortedRemoveDupsInt(n,ii);CHKERRQ(ierr); 1431db36a52SBarry Smith PetscFunctionReturn(0); 1441db36a52SBarry Smith } 1451db36a52SBarry Smith 14660e03357SMatthew G Knepley /*@ 147d9f1162dSJed Brown PetscFindInt - Finds integer in a sorted array of integers 14860e03357SMatthew G Knepley 14960e03357SMatthew G Knepley Not Collective 15060e03357SMatthew G Knepley 15160e03357SMatthew G Knepley Input Parameters: 15260e03357SMatthew G Knepley + key - the integer to locate 153d9f1162dSJed Brown . n - number of values in the array 15460e03357SMatthew G Knepley - ii - array of integers 15560e03357SMatthew G Knepley 15660e03357SMatthew G Knepley Output Parameter: 157d9f1162dSJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 15860e03357SMatthew G Knepley 15960e03357SMatthew G Knepley Level: intermediate 16060e03357SMatthew G Knepley 16160e03357SMatthew G Knepley Concepts: sorting^ints 16260e03357SMatthew G Knepley 16360e03357SMatthew G Knepley .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt() 16460e03357SMatthew G Knepley @*/ 165026ec6cbSMatthew G Knepley PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt ii[], PetscInt *loc) 16660e03357SMatthew G Knepley { 167d9f1162dSJed Brown PetscInt lo = 0,hi = n; 16860e03357SMatthew G Knepley 16960e03357SMatthew G Knepley PetscFunctionBegin; 17060e03357SMatthew G Knepley PetscValidPointer(loc,4); 17198781d82SMatthew G Knepley if (!n) {*loc = -1; PetscFunctionReturn(0);} 17298781d82SMatthew G Knepley PetscValidPointer(ii,3); 173d9f1162dSJed Brown while (hi - lo > 1) { 174d9f1162dSJed Brown PetscInt mid = lo + (hi - lo)/2; 175d9f1162dSJed Brown if (key < ii[mid]) hi = mid; 176d9f1162dSJed Brown else lo = mid; 17760e03357SMatthew G Knepley } 178d9f1162dSJed Brown *loc = key == ii[lo] ? lo : -(lo + (key > ii[lo]) + 1); 17960e03357SMatthew G Knepley PetscFunctionReturn(0); 18060e03357SMatthew G Knepley } 18160e03357SMatthew G Knepley 182d2aeb606SJed Brown /*@ 183d2aeb606SJed Brown PetscFindMPIInt - Finds MPI integer in a sorted array of integers 184d2aeb606SJed Brown 185d2aeb606SJed Brown Not Collective 186d2aeb606SJed Brown 187d2aeb606SJed Brown Input Parameters: 188d2aeb606SJed Brown + key - the integer to locate 189d2aeb606SJed Brown . n - number of values in the array 190d2aeb606SJed Brown - ii - array of integers 191d2aeb606SJed Brown 192d2aeb606SJed Brown Output Parameter: 193d2aeb606SJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 194d2aeb606SJed Brown 195d2aeb606SJed Brown Level: intermediate 196d2aeb606SJed Brown 197d2aeb606SJed Brown Concepts: sorting^ints 198d2aeb606SJed Brown 199d2aeb606SJed Brown .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt() 200d2aeb606SJed Brown @*/ 201d2aeb606SJed Brown PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt ii[], PetscInt *loc) 202d2aeb606SJed Brown { 203d2aeb606SJed Brown PetscInt lo = 0,hi = n; 204d2aeb606SJed Brown 205d2aeb606SJed Brown PetscFunctionBegin; 206d2aeb606SJed Brown PetscValidPointer(loc,4); 207d2aeb606SJed Brown if (!n) {*loc = -1; PetscFunctionReturn(0);} 208d2aeb606SJed Brown PetscValidPointer(ii,3); 209d2aeb606SJed Brown while (hi - lo > 1) { 210d2aeb606SJed Brown PetscInt mid = lo + (hi - lo)/2; 211d2aeb606SJed Brown if (key < ii[mid]) hi = mid; 212d2aeb606SJed Brown else lo = mid; 213d2aeb606SJed Brown } 214d2aeb606SJed Brown *loc = key == ii[lo] ? lo : -(lo + (key > ii[lo]) + 1); 215d2aeb606SJed Brown PetscFunctionReturn(0); 216d2aeb606SJed Brown } 2171db36a52SBarry Smith 218e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 219e5c89e4eSSatish Balay #define SWAP2(a,b,c,d,t) {t=a;a=b;b=t;t=c;c=d;d=t;} 220e5c89e4eSSatish Balay 221e5c89e4eSSatish Balay /* 222e5c89e4eSSatish Balay A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 223e5c89e4eSSatish Balay Assumes 0 origin for v, number of elements = right+1 (right is index of 224e5c89e4eSSatish Balay right-most member). 225e5c89e4eSSatish Balay */ 226e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithArray_Private(PetscInt *v,PetscInt *V,PetscInt right) 227e5c89e4eSSatish Balay { 228e5c89e4eSSatish Balay PetscErrorCode ierr; 229e5c89e4eSSatish Balay PetscInt i,vl,last,tmp; 230e5c89e4eSSatish Balay 231e5c89e4eSSatish Balay PetscFunctionBegin; 232e5c89e4eSSatish Balay if (right <= 1) { 233e5c89e4eSSatish Balay if (right == 1) { 234e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp); 235e5c89e4eSSatish Balay } 236e5c89e4eSSatish Balay PetscFunctionReturn(0); 237e5c89e4eSSatish Balay } 238e5c89e4eSSatish Balay SWAP2(v[0],v[right/2],V[0],V[right/2],tmp); 239e5c89e4eSSatish Balay vl = v[0]; 240e5c89e4eSSatish Balay last = 0; 241e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 242e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);} 243e5c89e4eSSatish Balay } 244e5c89e4eSSatish Balay SWAP2(v[0],v[last],V[0],V[last],tmp); 245e5c89e4eSSatish Balay ierr = PetscSortIntWithArray_Private(v,V,last-1);CHKERRQ(ierr); 246e5c89e4eSSatish Balay ierr = PetscSortIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 247e5c89e4eSSatish Balay PetscFunctionReturn(0); 248e5c89e4eSSatish Balay } 249e5c89e4eSSatish Balay 250e5c89e4eSSatish Balay /*@ 251e5c89e4eSSatish Balay PetscSortIntWithArray - Sorts an array of integers in place in increasing order; 252e5c89e4eSSatish Balay changes a second array to match the sorted first array. 253e5c89e4eSSatish Balay 254e5c89e4eSSatish Balay Not Collective 255e5c89e4eSSatish Balay 256e5c89e4eSSatish Balay Input Parameters: 257e5c89e4eSSatish Balay + n - number of values 258e5c89e4eSSatish Balay . i - array of integers 259e5c89e4eSSatish Balay - I - second array of integers 260e5c89e4eSSatish Balay 261e5c89e4eSSatish Balay Level: intermediate 262e5c89e4eSSatish Balay 263e5c89e4eSSatish Balay Concepts: sorting^ints with array 264e5c89e4eSSatish Balay 265e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 266e5c89e4eSSatish Balay @*/ 2677087cfbeSBarry Smith PetscErrorCode PetscSortIntWithArray(PetscInt n,PetscInt i[],PetscInt Ii[]) 268e5c89e4eSSatish Balay { 269e5c89e4eSSatish Balay PetscErrorCode ierr; 270e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 271e5c89e4eSSatish Balay 272e5c89e4eSSatish Balay PetscFunctionBegin; 273e5c89e4eSSatish Balay if (n<8) { 274e5c89e4eSSatish Balay for (k=0; k<n; k++) { 275e5c89e4eSSatish Balay ik = i[k]; 276e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 277e5c89e4eSSatish Balay if (ik > i[j]) { 278b7940d39SSatish Balay SWAP2(i[k],i[j],Ii[k],Ii[j],tmp); 279e5c89e4eSSatish Balay ik = i[k]; 280e5c89e4eSSatish Balay } 281e5c89e4eSSatish Balay } 282e5c89e4eSSatish Balay } 283e5c89e4eSSatish Balay } else { 284b7940d39SSatish Balay ierr = PetscSortIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr); 285e5c89e4eSSatish Balay } 286e5c89e4eSSatish Balay PetscFunctionReturn(0); 287e5c89e4eSSatish Balay } 288e5c89e4eSSatish Balay 289c1f0200aSDmitry Karpeev 290c1f0200aSDmitry Karpeev #define SWAP3(a,b,c,d,e,f,t) {t=a;a=b;b=t;t=c;c=d;d=t;t=e;e=f;f=t;} 291c1f0200aSDmitry Karpeev 292c1f0200aSDmitry Karpeev /* 293c1f0200aSDmitry Karpeev A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 294c1f0200aSDmitry Karpeev Assumes 0 origin for v, number of elements = right+1 (right is index of 295c1f0200aSDmitry Karpeev right-most member). 296c1f0200aSDmitry Karpeev */ 297d7aa01a8SHong Zhang static PetscErrorCode PetscSortIntWithArrayPair_Private(PetscInt *L,PetscInt *J, PetscInt *K, PetscInt right) 298c1f0200aSDmitry Karpeev { 299c1f0200aSDmitry Karpeev PetscErrorCode ierr; 300c1f0200aSDmitry Karpeev PetscInt i,vl,last,tmp; 301c1f0200aSDmitry Karpeev 302c1f0200aSDmitry Karpeev PetscFunctionBegin; 303c1f0200aSDmitry Karpeev if (right <= 1) { 304c1f0200aSDmitry Karpeev if (right == 1) { 305d7aa01a8SHong Zhang if (L[0] > L[1]) SWAP3(L[0],L[1],J[0],J[1],K[0],K[1],tmp); 306c1f0200aSDmitry Karpeev } 307c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 308c1f0200aSDmitry Karpeev } 309d7aa01a8SHong Zhang SWAP3(L[0],L[right/2],J[0],J[right/2],K[0],K[right/2],tmp); 310d7aa01a8SHong Zhang vl = L[0]; 311c1f0200aSDmitry Karpeev last = 0; 312c1f0200aSDmitry Karpeev for (i=1; i<=right; i++) { 313d7aa01a8SHong Zhang if (L[i] < vl) {last++; SWAP3(L[last],L[i],J[last],J[i],K[last],K[i],tmp);} 314c1f0200aSDmitry Karpeev } 315d7aa01a8SHong Zhang SWAP3(L[0],L[last],J[0],J[last],K[0],K[last],tmp); 316d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L,J,K,last-1);CHKERRQ(ierr); 317d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L+last+1,J+last+1,K+last+1,right-(last+1));CHKERRQ(ierr); 318c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 319c1f0200aSDmitry Karpeev } 320c1f0200aSDmitry Karpeev 321c1f0200aSDmitry Karpeev /*@ 322c1f0200aSDmitry Karpeev PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order; 323c1f0200aSDmitry Karpeev changes a pair of integer arrays to match the sorted first array. 324c1f0200aSDmitry Karpeev 325c1f0200aSDmitry Karpeev Not Collective 326c1f0200aSDmitry Karpeev 327c1f0200aSDmitry Karpeev Input Parameters: 328c1f0200aSDmitry Karpeev + n - number of values 329c1f0200aSDmitry Karpeev . I - array of integers 330c1f0200aSDmitry Karpeev . J - second array of integers (first array of the pair) 331c1f0200aSDmitry Karpeev - K - third array of integers (second array of the pair) 332c1f0200aSDmitry Karpeev 333c1f0200aSDmitry Karpeev Level: intermediate 334c1f0200aSDmitry Karpeev 335c1f0200aSDmitry Karpeev Concepts: sorting^ints with array pair 336c1f0200aSDmitry Karpeev 337c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray() 338c1f0200aSDmitry Karpeev @*/ 3396c2863d0SBarry Smith PetscErrorCode PetscSortIntWithArrayPair(PetscInt n,PetscInt L[],PetscInt J[], PetscInt K[]) 340c1f0200aSDmitry Karpeev { 341c1f0200aSDmitry Karpeev PetscErrorCode ierr; 342c1f0200aSDmitry Karpeev PetscInt j,k,tmp,ik; 343c1f0200aSDmitry Karpeev 344c1f0200aSDmitry Karpeev PetscFunctionBegin; 345c1f0200aSDmitry Karpeev if (n<8) { 346c1f0200aSDmitry Karpeev for (k=0; k<n; k++) { 347d7aa01a8SHong Zhang ik = L[k]; 348c1f0200aSDmitry Karpeev for (j=k+1; j<n; j++) { 349d7aa01a8SHong Zhang if (ik > L[j]) { 350d7aa01a8SHong Zhang SWAP3(L[k],L[j],J[k],J[j],K[k],K[j],tmp); 351d7aa01a8SHong Zhang ik = L[k]; 352c1f0200aSDmitry Karpeev } 353c1f0200aSDmitry Karpeev } 354c1f0200aSDmitry Karpeev } 355c1f0200aSDmitry Karpeev } else { 356d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L,J,K,n-1);CHKERRQ(ierr); 357c1f0200aSDmitry Karpeev } 358c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 359c1f0200aSDmitry Karpeev } 360c1f0200aSDmitry Karpeev 36117d7d925SStefano Zampini /* 36217d7d925SStefano Zampini A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 36317d7d925SStefano Zampini Assumes 0 origin for v, number of elements = right+1 (right is index of 36417d7d925SStefano Zampini right-most member). 36517d7d925SStefano Zampini */ 36617d7d925SStefano Zampini static void PetscSortMPIInt_Private(PetscMPIInt *v,PetscInt right) 36717d7d925SStefano Zampini { 36817d7d925SStefano Zampini PetscInt i,j; 36917d7d925SStefano Zampini PetscMPIInt pivot,tmp; 37017d7d925SStefano Zampini 37117d7d925SStefano Zampini if (right <= 1) { 37217d7d925SStefano Zampini if (right == 1) { 37317d7d925SStefano Zampini if (v[0] > v[1]) SWAP(v[0],v[1],tmp); 37417d7d925SStefano Zampini } 37517d7d925SStefano Zampini return; 37617d7d925SStefano Zampini } 37717d7d925SStefano Zampini i = MEDIAN(v,right); /* Choose a pivot */ 37817d7d925SStefano Zampini SWAP(v[0],v[i],tmp); /* Move it out of the way */ 37917d7d925SStefano Zampini pivot = v[0]; 38017d7d925SStefano Zampini for (i=0,j=right+1;; ) { 38117d7d925SStefano Zampini while (++i < j && v[i] <= pivot) ; /* Scan from the left */ 38217d7d925SStefano Zampini while (v[--j] > pivot) ; /* Scan from the right */ 38317d7d925SStefano Zampini if (i >= j) break; 38417d7d925SStefano Zampini SWAP(v[i],v[j],tmp); 38517d7d925SStefano Zampini } 38617d7d925SStefano Zampini SWAP(v[0],v[j],tmp); /* Put pivot back in place. */ 38717d7d925SStefano Zampini PetscSortMPIInt_Private(v,j-1); 38817d7d925SStefano Zampini PetscSortMPIInt_Private(v+j+1,right-(j+1)); 38917d7d925SStefano Zampini } 39017d7d925SStefano Zampini 39117d7d925SStefano Zampini /*@ 39217d7d925SStefano Zampini PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order. 39317d7d925SStefano Zampini 39417d7d925SStefano Zampini Not Collective 39517d7d925SStefano Zampini 39617d7d925SStefano Zampini Input Parameters: 39717d7d925SStefano Zampini + n - number of values 39817d7d925SStefano Zampini - i - array of integers 39917d7d925SStefano Zampini 40017d7d925SStefano Zampini Level: intermediate 40117d7d925SStefano Zampini 40217d7d925SStefano Zampini Concepts: sorting^ints 40317d7d925SStefano Zampini 40417d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation() 40517d7d925SStefano Zampini @*/ 40617d7d925SStefano Zampini PetscErrorCode PetscSortMPIInt(PetscInt n,PetscMPIInt i[]) 40717d7d925SStefano Zampini { 40817d7d925SStefano Zampini PetscInt j,k; 40917d7d925SStefano Zampini PetscMPIInt tmp,ik; 41017d7d925SStefano Zampini 41117d7d925SStefano Zampini PetscFunctionBegin; 41217d7d925SStefano Zampini if (n<8) { 41317d7d925SStefano Zampini for (k=0; k<n; k++) { 41417d7d925SStefano Zampini ik = i[k]; 41517d7d925SStefano Zampini for (j=k+1; j<n; j++) { 41617d7d925SStefano Zampini if (ik > i[j]) { 41717d7d925SStefano Zampini SWAP(i[k],i[j],tmp); 41817d7d925SStefano Zampini ik = i[k]; 41917d7d925SStefano Zampini } 42017d7d925SStefano Zampini } 42117d7d925SStefano Zampini } 422a297a907SKarl Rupp } else PetscSortMPIInt_Private(i,n-1); 42317d7d925SStefano Zampini PetscFunctionReturn(0); 42417d7d925SStefano Zampini } 42517d7d925SStefano Zampini 42617d7d925SStefano Zampini /*@ 42717d7d925SStefano Zampini PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries 42817d7d925SStefano Zampini 42917d7d925SStefano Zampini Not Collective 43017d7d925SStefano Zampini 43117d7d925SStefano Zampini Input Parameters: 43217d7d925SStefano Zampini + n - number of values 43317d7d925SStefano Zampini - ii - array of integers 43417d7d925SStefano Zampini 43517d7d925SStefano Zampini Output Parameter: 43617d7d925SStefano Zampini . n - number of non-redundant values 43717d7d925SStefano Zampini 43817d7d925SStefano Zampini Level: intermediate 43917d7d925SStefano Zampini 44017d7d925SStefano Zampini Concepts: sorting^ints 44117d7d925SStefano Zampini 44217d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt() 44317d7d925SStefano Zampini @*/ 44417d7d925SStefano Zampini PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt ii[]) 44517d7d925SStefano Zampini { 44617d7d925SStefano Zampini PetscErrorCode ierr; 44717d7d925SStefano Zampini PetscInt i,s = 0,N = *n, b = 0; 44817d7d925SStefano Zampini 44917d7d925SStefano Zampini PetscFunctionBegin; 45017d7d925SStefano Zampini ierr = PetscSortMPIInt(N,ii);CHKERRQ(ierr); 45117d7d925SStefano Zampini for (i=0; i<N-1; i++) { 452a297a907SKarl Rupp if (ii[b+s+1] != ii[b]) { 453a297a907SKarl Rupp ii[b+1] = ii[b+s+1]; b++; 454a297a907SKarl Rupp } else s++; 45517d7d925SStefano Zampini } 45617d7d925SStefano Zampini *n = N - s; 45717d7d925SStefano Zampini PetscFunctionReturn(0); 45817d7d925SStefano Zampini } 459c1f0200aSDmitry Karpeev 4604d615ea0SBarry Smith /* 4614d615ea0SBarry Smith A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 4624d615ea0SBarry Smith Assumes 0 origin for v, number of elements = right+1 (right is index of 4634d615ea0SBarry Smith right-most member). 4644d615ea0SBarry Smith */ 4654d615ea0SBarry Smith static PetscErrorCode PetscSortMPIIntWithArray_Private(PetscMPIInt *v,PetscMPIInt *V,PetscMPIInt right) 4664d615ea0SBarry Smith { 4674d615ea0SBarry Smith PetscErrorCode ierr; 4684d615ea0SBarry Smith PetscMPIInt i,vl,last,tmp; 4694d615ea0SBarry Smith 4704d615ea0SBarry Smith PetscFunctionBegin; 4714d615ea0SBarry Smith if (right <= 1) { 4724d615ea0SBarry Smith if (right == 1) { 4734d615ea0SBarry Smith if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp); 4744d615ea0SBarry Smith } 4754d615ea0SBarry Smith PetscFunctionReturn(0); 4764d615ea0SBarry Smith } 4774d615ea0SBarry Smith SWAP2(v[0],v[right/2],V[0],V[right/2],tmp); 4784d615ea0SBarry Smith vl = v[0]; 4794d615ea0SBarry Smith last = 0; 4804d615ea0SBarry Smith for (i=1; i<=right; i++) { 4814d615ea0SBarry Smith if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);} 4824d615ea0SBarry Smith } 4834d615ea0SBarry Smith SWAP2(v[0],v[last],V[0],V[last],tmp); 4844d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v,V,last-1);CHKERRQ(ierr); 4854d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 4864d615ea0SBarry Smith PetscFunctionReturn(0); 4874d615ea0SBarry Smith } 4884d615ea0SBarry Smith 4894d615ea0SBarry Smith /*@ 4904d615ea0SBarry Smith PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order; 4914d615ea0SBarry Smith changes a second array to match the sorted first array. 4924d615ea0SBarry Smith 4934d615ea0SBarry Smith Not Collective 4944d615ea0SBarry Smith 4954d615ea0SBarry Smith Input Parameters: 4964d615ea0SBarry Smith + n - number of values 4974d615ea0SBarry Smith . i - array of integers 4984d615ea0SBarry Smith - I - second array of integers 4994d615ea0SBarry Smith 5004d615ea0SBarry Smith Level: intermediate 5014d615ea0SBarry Smith 5024d615ea0SBarry Smith Concepts: sorting^ints with array 5034d615ea0SBarry Smith 5044d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 5054d615ea0SBarry Smith @*/ 5067087cfbeSBarry Smith PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt i[],PetscMPIInt Ii[]) 5074d615ea0SBarry Smith { 5084d615ea0SBarry Smith PetscErrorCode ierr; 5094d615ea0SBarry Smith PetscMPIInt j,k,tmp,ik; 5104d615ea0SBarry Smith 5114d615ea0SBarry Smith PetscFunctionBegin; 5124d615ea0SBarry Smith if (n<8) { 5134d615ea0SBarry Smith for (k=0; k<n; k++) { 5144d615ea0SBarry Smith ik = i[k]; 5154d615ea0SBarry Smith for (j=k+1; j<n; j++) { 5164d615ea0SBarry Smith if (ik > i[j]) { 5174d615ea0SBarry Smith SWAP2(i[k],i[j],Ii[k],Ii[j],tmp); 5184d615ea0SBarry Smith ik = i[k]; 5194d615ea0SBarry Smith } 5204d615ea0SBarry Smith } 5214d615ea0SBarry Smith } 5224d615ea0SBarry Smith } else { 5234d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr); 5244d615ea0SBarry Smith } 5254d615ea0SBarry Smith PetscFunctionReturn(0); 5264d615ea0SBarry Smith } 5274d615ea0SBarry Smith 528e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 529*5569a785SJunchao Zhang #define SWAP2MPIIntInt(a,b,c,d,t1,t2) {t1=a;a=b;b=t1;t2=c;c=d;d=t2;} 530*5569a785SJunchao Zhang 531*5569a785SJunchao Zhang /* 532*5569a785SJunchao Zhang A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 533*5569a785SJunchao Zhang Assumes 0 origin for v, number of elements = right+1 (right is index of 534*5569a785SJunchao Zhang right-most member). 535*5569a785SJunchao Zhang */ 536*5569a785SJunchao Zhang static PetscErrorCode PetscSortMPIIntWithIntArray_Private(PetscMPIInt *v,PetscInt *V,PetscMPIInt right) 537*5569a785SJunchao Zhang { 538*5569a785SJunchao Zhang PetscErrorCode ierr; 539*5569a785SJunchao Zhang PetscMPIInt i,vl,last,t1; 540*5569a785SJunchao Zhang PetscInt t2; 541*5569a785SJunchao Zhang 542*5569a785SJunchao Zhang PetscFunctionBegin; 543*5569a785SJunchao Zhang if (right <= 1) { 544*5569a785SJunchao Zhang if (right == 1) { 545*5569a785SJunchao Zhang if (v[0] > v[1]) SWAP2MPIIntInt(v[0],v[1],V[0],V[1],t1,t2); 546*5569a785SJunchao Zhang } 547*5569a785SJunchao Zhang PetscFunctionReturn(0); 548*5569a785SJunchao Zhang } 549*5569a785SJunchao Zhang SWAP2MPIIntInt(v[0],v[right/2],V[0],V[right/2],t1,t2); 550*5569a785SJunchao Zhang vl = v[0]; 551*5569a785SJunchao Zhang last = 0; 552*5569a785SJunchao Zhang for (i=1; i<=right; i++) { 553*5569a785SJunchao Zhang if (v[i] < vl) {last++; SWAP2MPIIntInt(v[last],v[i],V[last],V[i],t1,t2);} 554*5569a785SJunchao Zhang } 555*5569a785SJunchao Zhang SWAP2MPIIntInt(v[0],v[last],V[0],V[last],t1,t2); 556*5569a785SJunchao Zhang ierr = PetscSortMPIIntWithIntArray_Private(v,V,last-1);CHKERRQ(ierr); 557*5569a785SJunchao Zhang ierr = PetscSortMPIIntWithIntArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 558*5569a785SJunchao Zhang PetscFunctionReturn(0); 559*5569a785SJunchao Zhang } 560*5569a785SJunchao Zhang 561*5569a785SJunchao Zhang /*@ 562*5569a785SJunchao Zhang PetscSortMPIIntWithIntArray - Sorts an array of MPI integers in place in increasing order; 563*5569a785SJunchao Zhang changes a second array of Petsc intergers to match the sorted first array. 564*5569a785SJunchao Zhang 565*5569a785SJunchao Zhang Not Collective 566*5569a785SJunchao Zhang 567*5569a785SJunchao Zhang Input Parameters: 568*5569a785SJunchao Zhang + n - number of values 569*5569a785SJunchao Zhang . i - array of MPI integers 570*5569a785SJunchao Zhang - I - second array of Petsc integers 571*5569a785SJunchao Zhang 572*5569a785SJunchao Zhang Level: intermediate 573*5569a785SJunchao Zhang 574*5569a785SJunchao Zhang Notes: this routine is useful when one needs to sort MPI ranks with other integer arrays. 575*5569a785SJunchao Zhang 576*5569a785SJunchao Zhang Concepts: sorting^ints with array 577*5569a785SJunchao Zhang 578*5569a785SJunchao Zhang .seealso: PetscSortMPIIntWithArray() 579*5569a785SJunchao Zhang @*/ 580*5569a785SJunchao Zhang PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n,PetscMPIInt i[],PetscInt Ii[]) 581*5569a785SJunchao Zhang { 582*5569a785SJunchao Zhang PetscErrorCode ierr; 583*5569a785SJunchao Zhang PetscMPIInt j,k,t1,ik; 584*5569a785SJunchao Zhang PetscInt t2; 585*5569a785SJunchao Zhang 586*5569a785SJunchao Zhang PetscFunctionBegin; 587*5569a785SJunchao Zhang if (n<8) { 588*5569a785SJunchao Zhang for (k=0; k<n; k++) { 589*5569a785SJunchao Zhang ik = i[k]; 590*5569a785SJunchao Zhang for (j=k+1; j<n; j++) { 591*5569a785SJunchao Zhang if (ik > i[j]) { 592*5569a785SJunchao Zhang SWAP2MPIIntInt(i[k],i[j],Ii[k],Ii[j],t1,t2); 593*5569a785SJunchao Zhang ik = i[k]; 594*5569a785SJunchao Zhang } 595*5569a785SJunchao Zhang } 596*5569a785SJunchao Zhang } 597*5569a785SJunchao Zhang } else { 598*5569a785SJunchao Zhang ierr = PetscSortMPIIntWithIntArray_Private(i,Ii,n-1);CHKERRQ(ierr); 599*5569a785SJunchao Zhang } 600*5569a785SJunchao Zhang PetscFunctionReturn(0); 601*5569a785SJunchao Zhang } 602*5569a785SJunchao Zhang 603*5569a785SJunchao Zhang /* -----------------------------------------------------------------------*/ 604e5c89e4eSSatish Balay #define SWAP2IntScalar(a,b,c,d,t,ts) {t=a;a=b;b=t;ts=c;c=d;d=ts;} 605e5c89e4eSSatish Balay 606e5c89e4eSSatish Balay /* 607e5c89e4eSSatish Balay Modified from PetscSortIntWithArray_Private(). 608e5c89e4eSSatish Balay */ 609e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithScalarArray_Private(PetscInt *v,PetscScalar *V,PetscInt right) 610e5c89e4eSSatish Balay { 611e5c89e4eSSatish Balay PetscErrorCode ierr; 612e5c89e4eSSatish Balay PetscInt i,vl,last,tmp; 613e5c89e4eSSatish Balay PetscScalar stmp; 614e5c89e4eSSatish Balay 615e5c89e4eSSatish Balay PetscFunctionBegin; 616e5c89e4eSSatish Balay if (right <= 1) { 617e5c89e4eSSatish Balay if (right == 1) { 618e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP2IntScalar(v[0],v[1],V[0],V[1],tmp,stmp); 619e5c89e4eSSatish Balay } 620e5c89e4eSSatish Balay PetscFunctionReturn(0); 621e5c89e4eSSatish Balay } 622e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[right/2],V[0],V[right/2],tmp,stmp); 623e5c89e4eSSatish Balay vl = v[0]; 624e5c89e4eSSatish Balay last = 0; 625e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 626e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP2IntScalar(v[last],v[i],V[last],V[i],tmp,stmp);} 627e5c89e4eSSatish Balay } 628e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[last],V[0],V[last],tmp,stmp); 629e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v,V,last-1);CHKERRQ(ierr); 630e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 631e5c89e4eSSatish Balay PetscFunctionReturn(0); 632e5c89e4eSSatish Balay } 633e5c89e4eSSatish Balay 634e5c89e4eSSatish Balay /*@ 635e5c89e4eSSatish Balay PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order; 636e5c89e4eSSatish Balay changes a second SCALAR array to match the sorted first INTEGER array. 637e5c89e4eSSatish Balay 638e5c89e4eSSatish Balay Not Collective 639e5c89e4eSSatish Balay 640e5c89e4eSSatish Balay Input Parameters: 641e5c89e4eSSatish Balay + n - number of values 642e5c89e4eSSatish Balay . i - array of integers 643e5c89e4eSSatish Balay - I - second array of scalars 644e5c89e4eSSatish Balay 645e5c89e4eSSatish Balay Level: intermediate 646e5c89e4eSSatish Balay 647e5c89e4eSSatish Balay Concepts: sorting^ints with array 648e5c89e4eSSatish Balay 649e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 650e5c89e4eSSatish Balay @*/ 6517087cfbeSBarry Smith PetscErrorCode PetscSortIntWithScalarArray(PetscInt n,PetscInt i[],PetscScalar Ii[]) 652e5c89e4eSSatish Balay { 653e5c89e4eSSatish Balay PetscErrorCode ierr; 654e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 655e5c89e4eSSatish Balay PetscScalar stmp; 656e5c89e4eSSatish Balay 657e5c89e4eSSatish Balay PetscFunctionBegin; 658e5c89e4eSSatish Balay if (n<8) { 659e5c89e4eSSatish Balay for (k=0; k<n; k++) { 660e5c89e4eSSatish Balay ik = i[k]; 661e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 662e5c89e4eSSatish Balay if (ik > i[j]) { 663b7940d39SSatish Balay SWAP2IntScalar(i[k],i[j],Ii[k],Ii[j],tmp,stmp); 664e5c89e4eSSatish Balay ik = i[k]; 665e5c89e4eSSatish Balay } 666e5c89e4eSSatish Balay } 667e5c89e4eSSatish Balay } 668e5c89e4eSSatish Balay } else { 669b7940d39SSatish Balay ierr = PetscSortIntWithScalarArray_Private(i,Ii,n-1);CHKERRQ(ierr); 670e5c89e4eSSatish Balay } 671e5c89e4eSSatish Balay PetscFunctionReturn(0); 672e5c89e4eSSatish Balay } 673e5c89e4eSSatish Balay 674ece8f864SToby Isaac #define SWAP2IntData(a,b,c,d,t,td,siz) \ 675ece8f864SToby Isaac do { \ 676ece8f864SToby Isaac PetscErrorCode _ierr; \ 677ece8f864SToby Isaac t=a;a=b;b=t; \ 678ece8f864SToby Isaac _ierr = PetscMemcpy(td,c,siz);CHKERRQ(_ierr); \ 679ece8f864SToby Isaac _ierr = PetscMemcpy(c,d,siz);CHKERRQ(_ierr); \ 680ece8f864SToby Isaac _ierr = PetscMemcpy(d,td,siz);CHKERRQ(_ierr); \ 681ece8f864SToby Isaac } while(0) 68217df18f2SToby Isaac 68317df18f2SToby Isaac /* 68417df18f2SToby Isaac Modified from PetscSortIntWithArray_Private(). 68517df18f2SToby Isaac */ 68617df18f2SToby Isaac static PetscErrorCode PetscSortIntWithDataArray_Private(PetscInt *v,char *V,PetscInt right,size_t size,void *work) 68717df18f2SToby Isaac { 68817df18f2SToby Isaac PetscErrorCode ierr; 68917df18f2SToby Isaac PetscInt i,vl,last,tmp; 69017df18f2SToby Isaac 69117df18f2SToby Isaac PetscFunctionBegin; 69217df18f2SToby Isaac if (right <= 1) { 69317df18f2SToby Isaac if (right == 1) { 69417df18f2SToby Isaac if (v[0] > v[1]) SWAP2IntData(v[0],v[1],V,V+size,tmp,work,size); 69517df18f2SToby Isaac } 69617df18f2SToby Isaac PetscFunctionReturn(0); 69717df18f2SToby Isaac } 69817df18f2SToby Isaac SWAP2IntData(v[0],v[right/2],V,V+size*(right/2),tmp,work,size); 69917df18f2SToby Isaac vl = v[0]; 70017df18f2SToby Isaac last = 0; 70117df18f2SToby Isaac for (i=1; i<=right; i++) { 70217df18f2SToby Isaac if (v[i] < vl) {last++; SWAP2IntData(v[last],v[i],V+size*last,V+size*i,tmp,work,size);} 70317df18f2SToby Isaac } 70417df18f2SToby Isaac SWAP2IntData(v[0],v[last],V,V + size*last,tmp,work,size); 70517df18f2SToby Isaac ierr = PetscSortIntWithDataArray_Private(v,V,last-1,size,work);CHKERRQ(ierr); 70617df18f2SToby Isaac ierr = PetscSortIntWithDataArray_Private(v+last+1,V+size*(last+1),right-(last+1),size,work);CHKERRQ(ierr); 70717df18f2SToby Isaac PetscFunctionReturn(0); 70817df18f2SToby Isaac } 70917df18f2SToby Isaac 7106c2863d0SBarry Smith /*@C 71117df18f2SToby Isaac PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order; 71217df18f2SToby Isaac changes a second array to match the sorted first INTEGER array. Unlike other sort routines, the user must 71317df18f2SToby Isaac provide workspace (the size of an element in the data array) to use when sorting. 71417df18f2SToby Isaac 71517df18f2SToby Isaac Not Collective 71617df18f2SToby Isaac 71717df18f2SToby Isaac Input Parameters: 71817df18f2SToby Isaac + n - number of values 71917df18f2SToby Isaac . i - array of integers 72017df18f2SToby Isaac . Ii - second array of data 72117df18f2SToby Isaac . size - sizeof elements in the data array in bytes 72217df18f2SToby Isaac - work - workspace of "size" bytes used when sorting 72317df18f2SToby Isaac 72417df18f2SToby Isaac Level: intermediate 72517df18f2SToby Isaac 72617df18f2SToby Isaac Concepts: sorting^ints with array 72717df18f2SToby Isaac 72817df18f2SToby Isaac .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 72917df18f2SToby Isaac @*/ 73017df18f2SToby Isaac PetscErrorCode PetscSortIntWithDataArray(PetscInt n,PetscInt i[],void *Ii,size_t size,void *work) 73117df18f2SToby Isaac { 73217df18f2SToby Isaac char *V = (char *) Ii; 73317df18f2SToby Isaac PetscErrorCode ierr; 73417df18f2SToby Isaac PetscInt j,k,tmp,ik; 73517df18f2SToby Isaac 73617df18f2SToby Isaac PetscFunctionBegin; 73717df18f2SToby Isaac if (n<8) { 73817df18f2SToby Isaac for (k=0; k<n; k++) { 73917df18f2SToby Isaac ik = i[k]; 74017df18f2SToby Isaac for (j=k+1; j<n; j++) { 74117df18f2SToby Isaac if (ik > i[j]) { 74217df18f2SToby Isaac SWAP2IntData(i[k],i[j],V+size*k,V+size*j,tmp,work,size); 74317df18f2SToby Isaac ik = i[k]; 74417df18f2SToby Isaac } 74517df18f2SToby Isaac } 74617df18f2SToby Isaac } 74717df18f2SToby Isaac } else { 74817df18f2SToby Isaac ierr = PetscSortIntWithDataArray_Private(i,V,n-1,size,work);CHKERRQ(ierr); 74917df18f2SToby Isaac } 75017df18f2SToby Isaac PetscFunctionReturn(0); 75117df18f2SToby Isaac } 75217df18f2SToby Isaac 753b4301105SBarry Smith 75421e72a00SBarry Smith /*@ 75521e72a00SBarry Smith PetscMergeIntArray - Merges two SORTED integer arrays, removes duplicate elements. 75621e72a00SBarry Smith 75721e72a00SBarry Smith Not Collective 75821e72a00SBarry Smith 75921e72a00SBarry Smith Input Parameters: 76021e72a00SBarry Smith + an - number of values in the first array 76121e72a00SBarry Smith . aI - first sorted array of integers 76221e72a00SBarry Smith . bn - number of values in the second array 76321e72a00SBarry Smith - bI - second array of integers 76421e72a00SBarry Smith 76521e72a00SBarry Smith Output Parameters: 76621e72a00SBarry Smith + n - number of values in the merged array 7676c2863d0SBarry Smith - L - merged sorted array, this is allocated if an array is not provided 76821e72a00SBarry Smith 76921e72a00SBarry Smith Level: intermediate 77021e72a00SBarry Smith 77121e72a00SBarry Smith Concepts: merging^arrays 77221e72a00SBarry Smith 77321e72a00SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 77421e72a00SBarry Smith @*/ 7756c2863d0SBarry Smith PetscErrorCode PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L) 77621e72a00SBarry Smith { 77721e72a00SBarry Smith PetscErrorCode ierr; 77821e72a00SBarry Smith PetscInt *L_ = *L, ak, bk, k; 77921e72a00SBarry Smith 78021e72a00SBarry Smith if (!L_) { 78121e72a00SBarry Smith ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr); 78221e72a00SBarry Smith L_ = *L; 78321e72a00SBarry Smith } 78421e72a00SBarry Smith k = ak = bk = 0; 78521e72a00SBarry Smith while (ak < an && bk < bn) { 78621e72a00SBarry Smith if (aI[ak] == bI[bk]) { 78721e72a00SBarry Smith L_[k] = aI[ak]; 78821e72a00SBarry Smith ++ak; 78921e72a00SBarry Smith ++bk; 79021e72a00SBarry Smith ++k; 79121e72a00SBarry Smith } else if (aI[ak] < bI[bk]) { 79221e72a00SBarry Smith L_[k] = aI[ak]; 79321e72a00SBarry Smith ++ak; 79421e72a00SBarry Smith ++k; 79521e72a00SBarry Smith } else { 79621e72a00SBarry Smith L_[k] = bI[bk]; 79721e72a00SBarry Smith ++bk; 79821e72a00SBarry Smith ++k; 79921e72a00SBarry Smith } 80021e72a00SBarry Smith } 80121e72a00SBarry Smith if (ak < an) { 80221e72a00SBarry Smith ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr); 80321e72a00SBarry Smith k += (an-ak); 80421e72a00SBarry Smith } 80521e72a00SBarry Smith if (bk < bn) { 80621e72a00SBarry Smith ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr); 80721e72a00SBarry Smith k += (bn-bk); 80821e72a00SBarry Smith } 80921e72a00SBarry Smith *n = k; 81021e72a00SBarry Smith PetscFunctionReturn(0); 81121e72a00SBarry Smith } 812b4301105SBarry Smith 813c1f0200aSDmitry Karpeev /*@ 81421e72a00SBarry Smith PetscMergeIntArrayPair - Merges two SORTED integer arrays that share NO common values along with an additional array of integers. 815c1f0200aSDmitry Karpeev The additional arrays are the same length as sorted arrays and are merged 816c1f0200aSDmitry Karpeev in the order determined by the merging of the sorted pair. 817c1f0200aSDmitry Karpeev 818c1f0200aSDmitry Karpeev Not Collective 819c1f0200aSDmitry Karpeev 820c1f0200aSDmitry Karpeev Input Parameters: 821c1f0200aSDmitry Karpeev + an - number of values in the first array 822c1f0200aSDmitry Karpeev . aI - first sorted array of integers 823c1f0200aSDmitry Karpeev . aJ - first additional array of integers 824c1f0200aSDmitry Karpeev . bn - number of values in the second array 825c1f0200aSDmitry Karpeev . bI - second array of integers 826c1f0200aSDmitry Karpeev - bJ - second additional of integers 827c1f0200aSDmitry Karpeev 828c1f0200aSDmitry Karpeev Output Parameters: 829c1f0200aSDmitry Karpeev + n - number of values in the merged array (== an + bn) 83014c49006SJed Brown . L - merged sorted array 831c1f0200aSDmitry Karpeev - J - merged additional array 832c1f0200aSDmitry Karpeev 83395452b02SPatrick Sanan Notes: 83495452b02SPatrick 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 835c1f0200aSDmitry Karpeev Level: intermediate 836c1f0200aSDmitry Karpeev 837c1f0200aSDmitry Karpeev Concepts: merging^arrays 838c1f0200aSDmitry Karpeev 839c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 840c1f0200aSDmitry Karpeev @*/ 8416c2863d0SBarry 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) 842c1f0200aSDmitry Karpeev { 843c1f0200aSDmitry Karpeev PetscErrorCode ierr; 84470d8d27cSBarry Smith PetscInt n_, *L_, *J_, ak, bk, k; 845c1f0200aSDmitry Karpeev 84614c49006SJed Brown PetscFunctionBegin; 84770d8d27cSBarry Smith PetscValidIntPointer(L,8); 84870d8d27cSBarry Smith PetscValidIntPointer(J,9); 849c1f0200aSDmitry Karpeev n_ = an + bn; 850c1f0200aSDmitry Karpeev *n = n_; 85170d8d27cSBarry Smith if (!*L) { 852785e854fSJed Brown ierr = PetscMalloc1(n_, L);CHKERRQ(ierr); 85370d8d27cSBarry Smith } 854d7aa01a8SHong Zhang L_ = *L; 85570d8d27cSBarry Smith if (!*J) { 85670d8d27cSBarry Smith ierr = PetscMalloc1(n_, J);CHKERRQ(ierr); 857c1f0200aSDmitry Karpeev } 858c1f0200aSDmitry Karpeev J_ = *J; 859c1f0200aSDmitry Karpeev k = ak = bk = 0; 860c1f0200aSDmitry Karpeev while (ak < an && bk < bn) { 861c1f0200aSDmitry Karpeev if (aI[ak] <= bI[bk]) { 862d7aa01a8SHong Zhang L_[k] = aI[ak]; 863c1f0200aSDmitry Karpeev J_[k] = aJ[ak]; 864c1f0200aSDmitry Karpeev ++ak; 865c1f0200aSDmitry Karpeev ++k; 866a297a907SKarl Rupp } else { 867d7aa01a8SHong Zhang L_[k] = bI[bk]; 868c1f0200aSDmitry Karpeev J_[k] = bJ[bk]; 869c1f0200aSDmitry Karpeev ++bk; 870c1f0200aSDmitry Karpeev ++k; 871c1f0200aSDmitry Karpeev } 872c1f0200aSDmitry Karpeev } 873c1f0200aSDmitry Karpeev if (ak < an) { 874d7aa01a8SHong Zhang ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr); 875c1f0200aSDmitry Karpeev ierr = PetscMemcpy(J_+k,aJ+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr); 876c1f0200aSDmitry Karpeev k += (an-ak); 877c1f0200aSDmitry Karpeev } 878c1f0200aSDmitry Karpeev if (bk < bn) { 879d7aa01a8SHong Zhang ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr); 880c1f0200aSDmitry Karpeev ierr = PetscMemcpy(J_+k,bJ+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr); 881c1f0200aSDmitry Karpeev } 882c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 883c1f0200aSDmitry Karpeev } 884c1f0200aSDmitry Karpeev 885e498c390SJed Brown /*@ 886e498c390SJed Brown PetscMergeMPIIntArray - Merges two SORTED integer arrays. 887e498c390SJed Brown 888e498c390SJed Brown Not Collective 889e498c390SJed Brown 890e498c390SJed Brown Input Parameters: 891e498c390SJed Brown + an - number of values in the first array 892e498c390SJed Brown . aI - first sorted array of integers 893e498c390SJed Brown . bn - number of values in the second array 894e498c390SJed Brown - bI - second array of integers 895e498c390SJed Brown 896e498c390SJed Brown Output Parameters: 897e498c390SJed Brown + n - number of values in the merged array (<= an + bn) 898e498c390SJed Brown - L - merged sorted array, allocated if address of NULL pointer is passed 899e498c390SJed Brown 900e498c390SJed Brown Level: intermediate 901e498c390SJed Brown 902e498c390SJed Brown Concepts: merging^arrays 903e498c390SJed Brown 904e498c390SJed Brown .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 905e498c390SJed Brown @*/ 906e498c390SJed Brown PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L) 907e498c390SJed Brown { 908e498c390SJed Brown PetscErrorCode ierr; 909e498c390SJed Brown PetscInt ai,bi,k; 910e498c390SJed Brown 911e498c390SJed Brown PetscFunctionBegin; 912e498c390SJed Brown if (!*L) {ierr = PetscMalloc1((an+bn),L);CHKERRQ(ierr);} 913e498c390SJed Brown for (ai=0,bi=0,k=0; ai<an || bi<bn; ) { 914e498c390SJed Brown PetscInt t = -1; 915e498c390SJed Brown for ( ; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai]; 916e498c390SJed Brown for ( ; bi<bn && bI[bi] == t; bi++); 917e498c390SJed Brown for ( ; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi]; 918e498c390SJed Brown for ( ; ai<an && aI[ai] == t; ai++); 919e498c390SJed Brown } 920e498c390SJed Brown *n = k; 921e498c390SJed Brown PetscFunctionReturn(0); 922e498c390SJed Brown } 923e5c89e4eSSatish Balay 9246c2863d0SBarry Smith /*@C 92542eaa7f3SBarry Smith PetscProcessTree - Prepares tree data to be displayed graphically 92642eaa7f3SBarry Smith 92742eaa7f3SBarry Smith Not Collective 92842eaa7f3SBarry Smith 92942eaa7f3SBarry Smith Input Parameters: 93042eaa7f3SBarry Smith + n - number of values 93142eaa7f3SBarry Smith . mask - indicates those entries in the tree, location 0 is always masked 93242eaa7f3SBarry Smith - parentid - indicates the parent of each entry 93342eaa7f3SBarry Smith 93442eaa7f3SBarry Smith Output Parameters: 93542eaa7f3SBarry Smith + Nlevels - the number of levels 93642eaa7f3SBarry Smith . Level - for each node tells its level 93742eaa7f3SBarry Smith . Levelcnts - the number of nodes on each level 93842eaa7f3SBarry Smith . Idbylevel - a list of ids on each of the levels, first level followed by second etc 93942eaa7f3SBarry Smith - Column - for each id tells its column index 94042eaa7f3SBarry Smith 9416c2863d0SBarry Smith Level: developer 94242eaa7f3SBarry Smith 94395452b02SPatrick Sanan Notes: 94495452b02SPatrick Sanan This code is not currently used 94542eaa7f3SBarry Smith 94642eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation() 94742eaa7f3SBarry Smith @*/ 9487087cfbeSBarry Smith PetscErrorCode PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column) 94942eaa7f3SBarry Smith { 95042eaa7f3SBarry Smith PetscInt i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column; 95142eaa7f3SBarry Smith PetscErrorCode ierr; 952ace3abfcSBarry Smith PetscBool done = PETSC_FALSE; 95342eaa7f3SBarry Smith 95442eaa7f3SBarry Smith PetscFunctionBegin; 95542eaa7f3SBarry Smith if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set"); 95642eaa7f3SBarry Smith for (i=0; i<n; i++) { 95742eaa7f3SBarry Smith if (mask[i]) continue; 95842eaa7f3SBarry Smith if (parentid[i] == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent"); 95942eaa7f3SBarry Smith if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked"); 96042eaa7f3SBarry Smith } 96142eaa7f3SBarry Smith 96242eaa7f3SBarry Smith for (i=0; i<n; i++) { 96342eaa7f3SBarry Smith if (!mask[i]) nmask++; 96442eaa7f3SBarry Smith } 96542eaa7f3SBarry Smith 96642eaa7f3SBarry Smith /* determine the level in the tree of each node */ 9671795a4d1SJed Brown ierr = PetscCalloc1(n,&level);CHKERRQ(ierr); 968a297a907SKarl Rupp 96942eaa7f3SBarry Smith level[0] = 1; 97042eaa7f3SBarry Smith while (!done) { 97142eaa7f3SBarry Smith done = PETSC_TRUE; 97242eaa7f3SBarry Smith for (i=0; i<n; i++) { 97342eaa7f3SBarry Smith if (mask[i]) continue; 97442eaa7f3SBarry Smith if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1; 97542eaa7f3SBarry Smith else if (!level[i]) done = PETSC_FALSE; 97642eaa7f3SBarry Smith } 97742eaa7f3SBarry Smith } 97842eaa7f3SBarry Smith for (i=0; i<n; i++) { 97942eaa7f3SBarry Smith level[i]--; 98042eaa7f3SBarry Smith nlevels = PetscMax(nlevels,level[i]); 98142eaa7f3SBarry Smith } 98242eaa7f3SBarry Smith 98342eaa7f3SBarry Smith /* count the number of nodes on each level and its max */ 9841795a4d1SJed Brown ierr = PetscCalloc1(nlevels,&levelcnt);CHKERRQ(ierr); 98542eaa7f3SBarry Smith for (i=0; i<n; i++) { 98642eaa7f3SBarry Smith if (mask[i]) continue; 98742eaa7f3SBarry Smith levelcnt[level[i]-1]++; 98842eaa7f3SBarry Smith } 989a297a907SKarl Rupp for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]); 99042eaa7f3SBarry Smith 99142eaa7f3SBarry Smith /* for each level sort the ids by the parent id */ 992dcca6d9dSJed Brown ierr = PetscMalloc2(levelmax,&workid,levelmax,&workparentid);CHKERRQ(ierr); 993785e854fSJed Brown ierr = PetscMalloc1(nmask,&idbylevel);CHKERRQ(ierr); 99442eaa7f3SBarry Smith for (j=1; j<=nlevels;j++) { 99542eaa7f3SBarry Smith cnt = 0; 99642eaa7f3SBarry Smith for (i=0; i<n; i++) { 99742eaa7f3SBarry Smith if (mask[i]) continue; 99842eaa7f3SBarry Smith if (level[i] != j) continue; 99942eaa7f3SBarry Smith workid[cnt] = i; 100042eaa7f3SBarry Smith workparentid[cnt++] = parentid[i]; 100142eaa7f3SBarry Smith } 100242eaa7f3SBarry Smith /* PetscIntView(cnt,workparentid,0); 100342eaa7f3SBarry Smith PetscIntView(cnt,workid,0); 100442eaa7f3SBarry Smith ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr); 100542eaa7f3SBarry Smith PetscIntView(cnt,workparentid,0); 100642eaa7f3SBarry Smith PetscIntView(cnt,workid,0);*/ 100742eaa7f3SBarry Smith ierr = PetscMemcpy(idbylevel+tcnt,workid,cnt*sizeof(PetscInt));CHKERRQ(ierr); 100842eaa7f3SBarry Smith tcnt += cnt; 100942eaa7f3SBarry Smith } 101042eaa7f3SBarry Smith if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes"); 101142eaa7f3SBarry Smith ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr); 101242eaa7f3SBarry Smith 101342eaa7f3SBarry Smith /* for each node list its column */ 1014785e854fSJed Brown ierr = PetscMalloc1(n,&column);CHKERRQ(ierr); 101542eaa7f3SBarry Smith cnt = 0; 101642eaa7f3SBarry Smith for (j=0; j<nlevels; j++) { 101742eaa7f3SBarry Smith for (i=0; i<levelcnt[j]; i++) { 101842eaa7f3SBarry Smith column[idbylevel[cnt++]] = i; 101942eaa7f3SBarry Smith } 102042eaa7f3SBarry Smith } 102142eaa7f3SBarry Smith 102242eaa7f3SBarry Smith *Nlevels = nlevels; 102342eaa7f3SBarry Smith *Level = level; 102442eaa7f3SBarry Smith *Levelcnt = levelcnt; 102542eaa7f3SBarry Smith *Idbylevel = idbylevel; 102642eaa7f3SBarry Smith *Column = column; 102742eaa7f3SBarry Smith PetscFunctionReturn(0); 102842eaa7f3SBarry Smith } 1029