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 /*@ 86*22ab5688SLisandro Dalcin PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted input array 87*22ab5688SLisandro Dalcin 88*22ab5688SLisandro Dalcin Not Collective 89*22ab5688SLisandro Dalcin 90*22ab5688SLisandro Dalcin Input Parameters: 91*22ab5688SLisandro Dalcin + n - number of values 92*22ab5688SLisandro Dalcin - ii - sorted array of integers 93*22ab5688SLisandro Dalcin 94*22ab5688SLisandro Dalcin Output Parameter: 95*22ab5688SLisandro Dalcin . n - number of non-redundant values 96*22ab5688SLisandro Dalcin 97*22ab5688SLisandro Dalcin Level: intermediate 98*22ab5688SLisandro Dalcin 99*22ab5688SLisandro Dalcin Concepts: sorting^ints 100*22ab5688SLisandro Dalcin 101*22ab5688SLisandro Dalcin .seealso: PetscSortInt() 102*22ab5688SLisandro Dalcin @*/ 103*22ab5688SLisandro Dalcin PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n,PetscInt ii[]) 104*22ab5688SLisandro Dalcin { 105*22ab5688SLisandro Dalcin PetscInt i,s = 0,N = *n, b = 0; 106*22ab5688SLisandro Dalcin 107*22ab5688SLisandro Dalcin PetscFunctionBegin; 108*22ab5688SLisandro Dalcin for (i=0; i<N-1; i++) { 109*22ab5688SLisandro Dalcin if (PetscUnlikely(ii[b+s+1] < ii[b])) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Input array is not sorted"); 110*22ab5688SLisandro Dalcin if (ii[b+s+1] != ii[b]) { 111*22ab5688SLisandro Dalcin ii[b+1] = ii[b+s+1]; b++; 112*22ab5688SLisandro Dalcin } else s++; 113*22ab5688SLisandro Dalcin } 114*22ab5688SLisandro Dalcin *n = N - s; 115*22ab5688SLisandro Dalcin PetscFunctionReturn(0); 116*22ab5688SLisandro Dalcin } 117*22ab5688SLisandro Dalcin 118*22ab5688SLisandro 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 134*22ab5688SLisandro 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; 141*22ab5688SLisandro Dalcin ierr = PetscSortInt(*n,ii);CHKERRQ(ierr); 142*22ab5688SLisandro 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 #undef __FUNCT__ 183d2aeb606SJed Brown #define __FUNCT__ "PetscFindMPIInt" 184d2aeb606SJed Brown /*@ 185d2aeb606SJed Brown PetscFindMPIInt - Finds MPI integer in a sorted array of integers 186d2aeb606SJed Brown 187d2aeb606SJed Brown Not Collective 188d2aeb606SJed Brown 189d2aeb606SJed Brown Input Parameters: 190d2aeb606SJed Brown + key - the integer to locate 191d2aeb606SJed Brown . n - number of values in the array 192d2aeb606SJed Brown - ii - array of integers 193d2aeb606SJed Brown 194d2aeb606SJed Brown Output Parameter: 195d2aeb606SJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 196d2aeb606SJed Brown 197d2aeb606SJed Brown Level: intermediate 198d2aeb606SJed Brown 199d2aeb606SJed Brown Concepts: sorting^ints 200d2aeb606SJed Brown 201d2aeb606SJed Brown .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt() 202d2aeb606SJed Brown @*/ 203d2aeb606SJed Brown PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt ii[], PetscInt *loc) 204d2aeb606SJed Brown { 205d2aeb606SJed Brown PetscInt lo = 0,hi = n; 206d2aeb606SJed Brown 207d2aeb606SJed Brown PetscFunctionBegin; 208d2aeb606SJed Brown PetscValidPointer(loc,4); 209d2aeb606SJed Brown if (!n) {*loc = -1; PetscFunctionReturn(0);} 210d2aeb606SJed Brown PetscValidPointer(ii,3); 211d2aeb606SJed Brown while (hi - lo > 1) { 212d2aeb606SJed Brown PetscInt mid = lo + (hi - lo)/2; 213d2aeb606SJed Brown if (key < ii[mid]) hi = mid; 214d2aeb606SJed Brown else lo = mid; 215d2aeb606SJed Brown } 216d2aeb606SJed Brown *loc = key == ii[lo] ? lo : -(lo + (key > ii[lo]) + 1); 217d2aeb606SJed Brown PetscFunctionReturn(0); 218d2aeb606SJed Brown } 2191db36a52SBarry Smith 220e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 221e5c89e4eSSatish Balay #define SWAP2(a,b,c,d,t) {t=a;a=b;b=t;t=c;c=d;d=t;} 222e5c89e4eSSatish Balay 223e5c89e4eSSatish Balay /* 224e5c89e4eSSatish Balay A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 225e5c89e4eSSatish Balay Assumes 0 origin for v, number of elements = right+1 (right is index of 226e5c89e4eSSatish Balay right-most member). 227e5c89e4eSSatish Balay */ 228e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithArray_Private(PetscInt *v,PetscInt *V,PetscInt right) 229e5c89e4eSSatish Balay { 230e5c89e4eSSatish Balay PetscErrorCode ierr; 231e5c89e4eSSatish Balay PetscInt i,vl,last,tmp; 232e5c89e4eSSatish Balay 233e5c89e4eSSatish Balay PetscFunctionBegin; 234e5c89e4eSSatish Balay if (right <= 1) { 235e5c89e4eSSatish Balay if (right == 1) { 236e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp); 237e5c89e4eSSatish Balay } 238e5c89e4eSSatish Balay PetscFunctionReturn(0); 239e5c89e4eSSatish Balay } 240e5c89e4eSSatish Balay SWAP2(v[0],v[right/2],V[0],V[right/2],tmp); 241e5c89e4eSSatish Balay vl = v[0]; 242e5c89e4eSSatish Balay last = 0; 243e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 244e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);} 245e5c89e4eSSatish Balay } 246e5c89e4eSSatish Balay SWAP2(v[0],v[last],V[0],V[last],tmp); 247e5c89e4eSSatish Balay ierr = PetscSortIntWithArray_Private(v,V,last-1);CHKERRQ(ierr); 248e5c89e4eSSatish Balay ierr = PetscSortIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 249e5c89e4eSSatish Balay PetscFunctionReturn(0); 250e5c89e4eSSatish Balay } 251e5c89e4eSSatish Balay 252e5c89e4eSSatish Balay /*@ 253e5c89e4eSSatish Balay PetscSortIntWithArray - Sorts an array of integers in place in increasing order; 254e5c89e4eSSatish Balay changes a second array to match the sorted first array. 255e5c89e4eSSatish Balay 256e5c89e4eSSatish Balay Not Collective 257e5c89e4eSSatish Balay 258e5c89e4eSSatish Balay Input Parameters: 259e5c89e4eSSatish Balay + n - number of values 260e5c89e4eSSatish Balay . i - array of integers 261e5c89e4eSSatish Balay - I - second array of integers 262e5c89e4eSSatish Balay 263e5c89e4eSSatish Balay Level: intermediate 264e5c89e4eSSatish Balay 265e5c89e4eSSatish Balay Concepts: sorting^ints with array 266e5c89e4eSSatish Balay 267e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 268e5c89e4eSSatish Balay @*/ 2697087cfbeSBarry Smith PetscErrorCode PetscSortIntWithArray(PetscInt n,PetscInt i[],PetscInt Ii[]) 270e5c89e4eSSatish Balay { 271e5c89e4eSSatish Balay PetscErrorCode ierr; 272e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 273e5c89e4eSSatish Balay 274e5c89e4eSSatish Balay PetscFunctionBegin; 275e5c89e4eSSatish Balay if (n<8) { 276e5c89e4eSSatish Balay for (k=0; k<n; k++) { 277e5c89e4eSSatish Balay ik = i[k]; 278e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 279e5c89e4eSSatish Balay if (ik > i[j]) { 280b7940d39SSatish Balay SWAP2(i[k],i[j],Ii[k],Ii[j],tmp); 281e5c89e4eSSatish Balay ik = i[k]; 282e5c89e4eSSatish Balay } 283e5c89e4eSSatish Balay } 284e5c89e4eSSatish Balay } 285e5c89e4eSSatish Balay } else { 286b7940d39SSatish Balay ierr = PetscSortIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr); 287e5c89e4eSSatish Balay } 288e5c89e4eSSatish Balay PetscFunctionReturn(0); 289e5c89e4eSSatish Balay } 290e5c89e4eSSatish Balay 291c1f0200aSDmitry Karpeev 292c1f0200aSDmitry 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;} 293c1f0200aSDmitry Karpeev 294c1f0200aSDmitry Karpeev /* 295c1f0200aSDmitry Karpeev A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 296c1f0200aSDmitry Karpeev Assumes 0 origin for v, number of elements = right+1 (right is index of 297c1f0200aSDmitry Karpeev right-most member). 298c1f0200aSDmitry Karpeev */ 299d7aa01a8SHong Zhang static PetscErrorCode PetscSortIntWithArrayPair_Private(PetscInt *L,PetscInt *J, PetscInt *K, PetscInt right) 300c1f0200aSDmitry Karpeev { 301c1f0200aSDmitry Karpeev PetscErrorCode ierr; 302c1f0200aSDmitry Karpeev PetscInt i,vl,last,tmp; 303c1f0200aSDmitry Karpeev 304c1f0200aSDmitry Karpeev PetscFunctionBegin; 305c1f0200aSDmitry Karpeev if (right <= 1) { 306c1f0200aSDmitry Karpeev if (right == 1) { 307d7aa01a8SHong Zhang if (L[0] > L[1]) SWAP3(L[0],L[1],J[0],J[1],K[0],K[1],tmp); 308c1f0200aSDmitry Karpeev } 309c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 310c1f0200aSDmitry Karpeev } 311d7aa01a8SHong Zhang SWAP3(L[0],L[right/2],J[0],J[right/2],K[0],K[right/2],tmp); 312d7aa01a8SHong Zhang vl = L[0]; 313c1f0200aSDmitry Karpeev last = 0; 314c1f0200aSDmitry Karpeev for (i=1; i<=right; i++) { 315d7aa01a8SHong Zhang if (L[i] < vl) {last++; SWAP3(L[last],L[i],J[last],J[i],K[last],K[i],tmp);} 316c1f0200aSDmitry Karpeev } 317d7aa01a8SHong Zhang SWAP3(L[0],L[last],J[0],J[last],K[0],K[last],tmp); 318d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L,J,K,last-1);CHKERRQ(ierr); 319d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L+last+1,J+last+1,K+last+1,right-(last+1));CHKERRQ(ierr); 320c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 321c1f0200aSDmitry Karpeev } 322c1f0200aSDmitry Karpeev 323c1f0200aSDmitry Karpeev /*@ 324c1f0200aSDmitry Karpeev PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order; 325c1f0200aSDmitry Karpeev changes a pair of integer arrays to match the sorted first array. 326c1f0200aSDmitry Karpeev 327c1f0200aSDmitry Karpeev Not Collective 328c1f0200aSDmitry Karpeev 329c1f0200aSDmitry Karpeev Input Parameters: 330c1f0200aSDmitry Karpeev + n - number of values 331c1f0200aSDmitry Karpeev . I - array of integers 332c1f0200aSDmitry Karpeev . J - second array of integers (first array of the pair) 333c1f0200aSDmitry Karpeev - K - third array of integers (second array of the pair) 334c1f0200aSDmitry Karpeev 335c1f0200aSDmitry Karpeev Level: intermediate 336c1f0200aSDmitry Karpeev 337c1f0200aSDmitry Karpeev Concepts: sorting^ints with array pair 338c1f0200aSDmitry Karpeev 339c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray() 340c1f0200aSDmitry Karpeev @*/ 3416c2863d0SBarry Smith PetscErrorCode PetscSortIntWithArrayPair(PetscInt n,PetscInt L[],PetscInt J[], PetscInt K[]) 342c1f0200aSDmitry Karpeev { 343c1f0200aSDmitry Karpeev PetscErrorCode ierr; 344c1f0200aSDmitry Karpeev PetscInt j,k,tmp,ik; 345c1f0200aSDmitry Karpeev 346c1f0200aSDmitry Karpeev PetscFunctionBegin; 347c1f0200aSDmitry Karpeev if (n<8) { 348c1f0200aSDmitry Karpeev for (k=0; k<n; k++) { 349d7aa01a8SHong Zhang ik = L[k]; 350c1f0200aSDmitry Karpeev for (j=k+1; j<n; j++) { 351d7aa01a8SHong Zhang if (ik > L[j]) { 352d7aa01a8SHong Zhang SWAP3(L[k],L[j],J[k],J[j],K[k],K[j],tmp); 353d7aa01a8SHong Zhang ik = L[k]; 354c1f0200aSDmitry Karpeev } 355c1f0200aSDmitry Karpeev } 356c1f0200aSDmitry Karpeev } 357c1f0200aSDmitry Karpeev } else { 358d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L,J,K,n-1);CHKERRQ(ierr); 359c1f0200aSDmitry Karpeev } 360c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 361c1f0200aSDmitry Karpeev } 362c1f0200aSDmitry Karpeev 36317d7d925SStefano Zampini /* 36417d7d925SStefano Zampini A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 36517d7d925SStefano Zampini Assumes 0 origin for v, number of elements = right+1 (right is index of 36617d7d925SStefano Zampini right-most member). 36717d7d925SStefano Zampini */ 36817d7d925SStefano Zampini static void PetscSortMPIInt_Private(PetscMPIInt *v,PetscInt right) 36917d7d925SStefano Zampini { 37017d7d925SStefano Zampini PetscInt i,j; 37117d7d925SStefano Zampini PetscMPIInt pivot,tmp; 37217d7d925SStefano Zampini 37317d7d925SStefano Zampini if (right <= 1) { 37417d7d925SStefano Zampini if (right == 1) { 37517d7d925SStefano Zampini if (v[0] > v[1]) SWAP(v[0],v[1],tmp); 37617d7d925SStefano Zampini } 37717d7d925SStefano Zampini return; 37817d7d925SStefano Zampini } 37917d7d925SStefano Zampini i = MEDIAN(v,right); /* Choose a pivot */ 38017d7d925SStefano Zampini SWAP(v[0],v[i],tmp); /* Move it out of the way */ 38117d7d925SStefano Zampini pivot = v[0]; 38217d7d925SStefano Zampini for (i=0,j=right+1;; ) { 38317d7d925SStefano Zampini while (++i < j && v[i] <= pivot) ; /* Scan from the left */ 38417d7d925SStefano Zampini while (v[--j] > pivot) ; /* Scan from the right */ 38517d7d925SStefano Zampini if (i >= j) break; 38617d7d925SStefano Zampini SWAP(v[i],v[j],tmp); 38717d7d925SStefano Zampini } 38817d7d925SStefano Zampini SWAP(v[0],v[j],tmp); /* Put pivot back in place. */ 38917d7d925SStefano Zampini PetscSortMPIInt_Private(v,j-1); 39017d7d925SStefano Zampini PetscSortMPIInt_Private(v+j+1,right-(j+1)); 39117d7d925SStefano Zampini } 39217d7d925SStefano Zampini 39317d7d925SStefano Zampini /*@ 39417d7d925SStefano Zampini PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order. 39517d7d925SStefano Zampini 39617d7d925SStefano Zampini Not Collective 39717d7d925SStefano Zampini 39817d7d925SStefano Zampini Input Parameters: 39917d7d925SStefano Zampini + n - number of values 40017d7d925SStefano Zampini - i - array of integers 40117d7d925SStefano Zampini 40217d7d925SStefano Zampini Level: intermediate 40317d7d925SStefano Zampini 40417d7d925SStefano Zampini Concepts: sorting^ints 40517d7d925SStefano Zampini 40617d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation() 40717d7d925SStefano Zampini @*/ 40817d7d925SStefano Zampini PetscErrorCode PetscSortMPIInt(PetscInt n,PetscMPIInt i[]) 40917d7d925SStefano Zampini { 41017d7d925SStefano Zampini PetscInt j,k; 41117d7d925SStefano Zampini PetscMPIInt tmp,ik; 41217d7d925SStefano Zampini 41317d7d925SStefano Zampini PetscFunctionBegin; 41417d7d925SStefano Zampini if (n<8) { 41517d7d925SStefano Zampini for (k=0; k<n; k++) { 41617d7d925SStefano Zampini ik = i[k]; 41717d7d925SStefano Zampini for (j=k+1; j<n; j++) { 41817d7d925SStefano Zampini if (ik > i[j]) { 41917d7d925SStefano Zampini SWAP(i[k],i[j],tmp); 42017d7d925SStefano Zampini ik = i[k]; 42117d7d925SStefano Zampini } 42217d7d925SStefano Zampini } 42317d7d925SStefano Zampini } 424a297a907SKarl Rupp } else PetscSortMPIInt_Private(i,n-1); 42517d7d925SStefano Zampini PetscFunctionReturn(0); 42617d7d925SStefano Zampini } 42717d7d925SStefano Zampini 42817d7d925SStefano Zampini /*@ 42917d7d925SStefano Zampini PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries 43017d7d925SStefano Zampini 43117d7d925SStefano Zampini Not Collective 43217d7d925SStefano Zampini 43317d7d925SStefano Zampini Input Parameters: 43417d7d925SStefano Zampini + n - number of values 43517d7d925SStefano Zampini - ii - array of integers 43617d7d925SStefano Zampini 43717d7d925SStefano Zampini Output Parameter: 43817d7d925SStefano Zampini . n - number of non-redundant values 43917d7d925SStefano Zampini 44017d7d925SStefano Zampini Level: intermediate 44117d7d925SStefano Zampini 44217d7d925SStefano Zampini Concepts: sorting^ints 44317d7d925SStefano Zampini 44417d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt() 44517d7d925SStefano Zampini @*/ 44617d7d925SStefano Zampini PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt ii[]) 44717d7d925SStefano Zampini { 44817d7d925SStefano Zampini PetscErrorCode ierr; 44917d7d925SStefano Zampini PetscInt i,s = 0,N = *n, b = 0; 45017d7d925SStefano Zampini 45117d7d925SStefano Zampini PetscFunctionBegin; 45217d7d925SStefano Zampini ierr = PetscSortMPIInt(N,ii);CHKERRQ(ierr); 45317d7d925SStefano Zampini for (i=0; i<N-1; i++) { 454a297a907SKarl Rupp if (ii[b+s+1] != ii[b]) { 455a297a907SKarl Rupp ii[b+1] = ii[b+s+1]; b++; 456a297a907SKarl Rupp } else s++; 45717d7d925SStefano Zampini } 45817d7d925SStefano Zampini *n = N - s; 45917d7d925SStefano Zampini PetscFunctionReturn(0); 46017d7d925SStefano Zampini } 461c1f0200aSDmitry Karpeev 4624d615ea0SBarry Smith /* 4634d615ea0SBarry Smith A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 4644d615ea0SBarry Smith Assumes 0 origin for v, number of elements = right+1 (right is index of 4654d615ea0SBarry Smith right-most member). 4664d615ea0SBarry Smith */ 4674d615ea0SBarry Smith static PetscErrorCode PetscSortMPIIntWithArray_Private(PetscMPIInt *v,PetscMPIInt *V,PetscMPIInt right) 4684d615ea0SBarry Smith { 4694d615ea0SBarry Smith PetscErrorCode ierr; 4704d615ea0SBarry Smith PetscMPIInt i,vl,last,tmp; 4714d615ea0SBarry Smith 4724d615ea0SBarry Smith PetscFunctionBegin; 4734d615ea0SBarry Smith if (right <= 1) { 4744d615ea0SBarry Smith if (right == 1) { 4754d615ea0SBarry Smith if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp); 4764d615ea0SBarry Smith } 4774d615ea0SBarry Smith PetscFunctionReturn(0); 4784d615ea0SBarry Smith } 4794d615ea0SBarry Smith SWAP2(v[0],v[right/2],V[0],V[right/2],tmp); 4804d615ea0SBarry Smith vl = v[0]; 4814d615ea0SBarry Smith last = 0; 4824d615ea0SBarry Smith for (i=1; i<=right; i++) { 4834d615ea0SBarry Smith if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);} 4844d615ea0SBarry Smith } 4854d615ea0SBarry Smith SWAP2(v[0],v[last],V[0],V[last],tmp); 4864d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v,V,last-1);CHKERRQ(ierr); 4874d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 4884d615ea0SBarry Smith PetscFunctionReturn(0); 4894d615ea0SBarry Smith } 4904d615ea0SBarry Smith 4914d615ea0SBarry Smith /*@ 4924d615ea0SBarry Smith PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order; 4934d615ea0SBarry Smith changes a second array to match the sorted first array. 4944d615ea0SBarry Smith 4954d615ea0SBarry Smith Not Collective 4964d615ea0SBarry Smith 4974d615ea0SBarry Smith Input Parameters: 4984d615ea0SBarry Smith + n - number of values 4994d615ea0SBarry Smith . i - array of integers 5004d615ea0SBarry Smith - I - second array of integers 5014d615ea0SBarry Smith 5024d615ea0SBarry Smith Level: intermediate 5034d615ea0SBarry Smith 5044d615ea0SBarry Smith Concepts: sorting^ints with array 5054d615ea0SBarry Smith 5064d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 5074d615ea0SBarry Smith @*/ 5087087cfbeSBarry Smith PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt i[],PetscMPIInt Ii[]) 5094d615ea0SBarry Smith { 5104d615ea0SBarry Smith PetscErrorCode ierr; 5114d615ea0SBarry Smith PetscMPIInt j,k,tmp,ik; 5124d615ea0SBarry Smith 5134d615ea0SBarry Smith PetscFunctionBegin; 5144d615ea0SBarry Smith if (n<8) { 5154d615ea0SBarry Smith for (k=0; k<n; k++) { 5164d615ea0SBarry Smith ik = i[k]; 5174d615ea0SBarry Smith for (j=k+1; j<n; j++) { 5184d615ea0SBarry Smith if (ik > i[j]) { 5194d615ea0SBarry Smith SWAP2(i[k],i[j],Ii[k],Ii[j],tmp); 5204d615ea0SBarry Smith ik = i[k]; 5214d615ea0SBarry Smith } 5224d615ea0SBarry Smith } 5234d615ea0SBarry Smith } 5244d615ea0SBarry Smith } else { 5254d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr); 5264d615ea0SBarry Smith } 5274d615ea0SBarry Smith PetscFunctionReturn(0); 5284d615ea0SBarry Smith } 5294d615ea0SBarry Smith 530e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 531e5c89e4eSSatish Balay #define SWAP2IntScalar(a,b,c,d,t,ts) {t=a;a=b;b=t;ts=c;c=d;d=ts;} 532e5c89e4eSSatish Balay 533e5c89e4eSSatish Balay /* 534e5c89e4eSSatish Balay Modified from PetscSortIntWithArray_Private(). 535e5c89e4eSSatish Balay */ 536e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithScalarArray_Private(PetscInt *v,PetscScalar *V,PetscInt right) 537e5c89e4eSSatish Balay { 538e5c89e4eSSatish Balay PetscErrorCode ierr; 539e5c89e4eSSatish Balay PetscInt i,vl,last,tmp; 540e5c89e4eSSatish Balay PetscScalar stmp; 541e5c89e4eSSatish Balay 542e5c89e4eSSatish Balay PetscFunctionBegin; 543e5c89e4eSSatish Balay if (right <= 1) { 544e5c89e4eSSatish Balay if (right == 1) { 545e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP2IntScalar(v[0],v[1],V[0],V[1],tmp,stmp); 546e5c89e4eSSatish Balay } 547e5c89e4eSSatish Balay PetscFunctionReturn(0); 548e5c89e4eSSatish Balay } 549e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[right/2],V[0],V[right/2],tmp,stmp); 550e5c89e4eSSatish Balay vl = v[0]; 551e5c89e4eSSatish Balay last = 0; 552e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 553e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP2IntScalar(v[last],v[i],V[last],V[i],tmp,stmp);} 554e5c89e4eSSatish Balay } 555e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[last],V[0],V[last],tmp,stmp); 556e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v,V,last-1);CHKERRQ(ierr); 557e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 558e5c89e4eSSatish Balay PetscFunctionReturn(0); 559e5c89e4eSSatish Balay } 560e5c89e4eSSatish Balay 561e5c89e4eSSatish Balay /*@ 562e5c89e4eSSatish Balay PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order; 563e5c89e4eSSatish Balay changes a second SCALAR array to match the sorted first INTEGER array. 564e5c89e4eSSatish Balay 565e5c89e4eSSatish Balay Not Collective 566e5c89e4eSSatish Balay 567e5c89e4eSSatish Balay Input Parameters: 568e5c89e4eSSatish Balay + n - number of values 569e5c89e4eSSatish Balay . i - array of integers 570e5c89e4eSSatish Balay - I - second array of scalars 571e5c89e4eSSatish Balay 572e5c89e4eSSatish Balay Level: intermediate 573e5c89e4eSSatish Balay 574e5c89e4eSSatish Balay Concepts: sorting^ints with array 575e5c89e4eSSatish Balay 576e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 577e5c89e4eSSatish Balay @*/ 5787087cfbeSBarry Smith PetscErrorCode PetscSortIntWithScalarArray(PetscInt n,PetscInt i[],PetscScalar Ii[]) 579e5c89e4eSSatish Balay { 580e5c89e4eSSatish Balay PetscErrorCode ierr; 581e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 582e5c89e4eSSatish Balay PetscScalar stmp; 583e5c89e4eSSatish Balay 584e5c89e4eSSatish Balay PetscFunctionBegin; 585e5c89e4eSSatish Balay if (n<8) { 586e5c89e4eSSatish Balay for (k=0; k<n; k++) { 587e5c89e4eSSatish Balay ik = i[k]; 588e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 589e5c89e4eSSatish Balay if (ik > i[j]) { 590b7940d39SSatish Balay SWAP2IntScalar(i[k],i[j],Ii[k],Ii[j],tmp,stmp); 591e5c89e4eSSatish Balay ik = i[k]; 592e5c89e4eSSatish Balay } 593e5c89e4eSSatish Balay } 594e5c89e4eSSatish Balay } 595e5c89e4eSSatish Balay } else { 596b7940d39SSatish Balay ierr = PetscSortIntWithScalarArray_Private(i,Ii,n-1);CHKERRQ(ierr); 597e5c89e4eSSatish Balay } 598e5c89e4eSSatish Balay PetscFunctionReturn(0); 599e5c89e4eSSatish Balay } 600e5c89e4eSSatish Balay 601ece8f864SToby Isaac #define SWAP2IntData(a,b,c,d,t,td,siz) \ 602ece8f864SToby Isaac do { \ 603ece8f864SToby Isaac PetscErrorCode _ierr; \ 604ece8f864SToby Isaac t=a;a=b;b=t; \ 605ece8f864SToby Isaac _ierr = PetscMemcpy(td,c,siz);CHKERRQ(_ierr); \ 606ece8f864SToby Isaac _ierr = PetscMemcpy(c,d,siz);CHKERRQ(_ierr); \ 607ece8f864SToby Isaac _ierr = PetscMemcpy(d,td,siz);CHKERRQ(_ierr); \ 608ece8f864SToby Isaac } while(0) 60917df18f2SToby Isaac 61017df18f2SToby Isaac /* 61117df18f2SToby Isaac Modified from PetscSortIntWithArray_Private(). 61217df18f2SToby Isaac */ 61317df18f2SToby Isaac static PetscErrorCode PetscSortIntWithDataArray_Private(PetscInt *v,char *V,PetscInt right,size_t size,void *work) 61417df18f2SToby Isaac { 61517df18f2SToby Isaac PetscErrorCode ierr; 61617df18f2SToby Isaac PetscInt i,vl,last,tmp; 61717df18f2SToby Isaac 61817df18f2SToby Isaac PetscFunctionBegin; 61917df18f2SToby Isaac if (right <= 1) { 62017df18f2SToby Isaac if (right == 1) { 62117df18f2SToby Isaac if (v[0] > v[1]) SWAP2IntData(v[0],v[1],V,V+size,tmp,work,size); 62217df18f2SToby Isaac } 62317df18f2SToby Isaac PetscFunctionReturn(0); 62417df18f2SToby Isaac } 62517df18f2SToby Isaac SWAP2IntData(v[0],v[right/2],V,V+size*(right/2),tmp,work,size); 62617df18f2SToby Isaac vl = v[0]; 62717df18f2SToby Isaac last = 0; 62817df18f2SToby Isaac for (i=1; i<=right; i++) { 62917df18f2SToby Isaac if (v[i] < vl) {last++; SWAP2IntData(v[last],v[i],V+size*last,V+size*i,tmp,work,size);} 63017df18f2SToby Isaac } 63117df18f2SToby Isaac SWAP2IntData(v[0],v[last],V,V + size*last,tmp,work,size); 63217df18f2SToby Isaac ierr = PetscSortIntWithDataArray_Private(v,V,last-1,size,work);CHKERRQ(ierr); 63317df18f2SToby Isaac ierr = PetscSortIntWithDataArray_Private(v+last+1,V+size*(last+1),right-(last+1),size,work);CHKERRQ(ierr); 63417df18f2SToby Isaac PetscFunctionReturn(0); 63517df18f2SToby Isaac } 63617df18f2SToby Isaac 6376c2863d0SBarry Smith /*@C 63817df18f2SToby Isaac PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order; 63917df18f2SToby Isaac changes a second array to match the sorted first INTEGER array. Unlike other sort routines, the user must 64017df18f2SToby Isaac provide workspace (the size of an element in the data array) to use when sorting. 64117df18f2SToby Isaac 64217df18f2SToby Isaac Not Collective 64317df18f2SToby Isaac 64417df18f2SToby Isaac Input Parameters: 64517df18f2SToby Isaac + n - number of values 64617df18f2SToby Isaac . i - array of integers 64717df18f2SToby Isaac . Ii - second array of data 64817df18f2SToby Isaac . size - sizeof elements in the data array in bytes 64917df18f2SToby Isaac - work - workspace of "size" bytes used when sorting 65017df18f2SToby Isaac 65117df18f2SToby Isaac Level: intermediate 65217df18f2SToby Isaac 65317df18f2SToby Isaac Concepts: sorting^ints with array 65417df18f2SToby Isaac 65517df18f2SToby Isaac .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 65617df18f2SToby Isaac @*/ 65717df18f2SToby Isaac PetscErrorCode PetscSortIntWithDataArray(PetscInt n,PetscInt i[],void *Ii,size_t size,void *work) 65817df18f2SToby Isaac { 65917df18f2SToby Isaac char *V = (char *) Ii; 66017df18f2SToby Isaac PetscErrorCode ierr; 66117df18f2SToby Isaac PetscInt j,k,tmp,ik; 66217df18f2SToby Isaac 66317df18f2SToby Isaac PetscFunctionBegin; 66417df18f2SToby Isaac if (n<8) { 66517df18f2SToby Isaac for (k=0; k<n; k++) { 66617df18f2SToby Isaac ik = i[k]; 66717df18f2SToby Isaac for (j=k+1; j<n; j++) { 66817df18f2SToby Isaac if (ik > i[j]) { 66917df18f2SToby Isaac SWAP2IntData(i[k],i[j],V+size*k,V+size*j,tmp,work,size); 67017df18f2SToby Isaac ik = i[k]; 67117df18f2SToby Isaac } 67217df18f2SToby Isaac } 67317df18f2SToby Isaac } 67417df18f2SToby Isaac } else { 67517df18f2SToby Isaac ierr = PetscSortIntWithDataArray_Private(i,V,n-1,size,work);CHKERRQ(ierr); 67617df18f2SToby Isaac } 67717df18f2SToby Isaac PetscFunctionReturn(0); 67817df18f2SToby Isaac } 67917df18f2SToby Isaac 680b4301105SBarry Smith 68121e72a00SBarry Smith /*@ 68221e72a00SBarry Smith PetscMergeIntArray - Merges two SORTED integer arrays, removes duplicate elements. 68321e72a00SBarry Smith 68421e72a00SBarry Smith Not Collective 68521e72a00SBarry Smith 68621e72a00SBarry Smith Input Parameters: 68721e72a00SBarry Smith + an - number of values in the first array 68821e72a00SBarry Smith . aI - first sorted array of integers 68921e72a00SBarry Smith . bn - number of values in the second array 69021e72a00SBarry Smith - bI - second array of integers 69121e72a00SBarry Smith 69221e72a00SBarry Smith Output Parameters: 69321e72a00SBarry Smith + n - number of values in the merged array 6946c2863d0SBarry Smith - L - merged sorted array, this is allocated if an array is not provided 69521e72a00SBarry Smith 69621e72a00SBarry Smith Level: intermediate 69721e72a00SBarry Smith 69821e72a00SBarry Smith Concepts: merging^arrays 69921e72a00SBarry Smith 70021e72a00SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 70121e72a00SBarry Smith @*/ 7026c2863d0SBarry Smith PetscErrorCode PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L) 70321e72a00SBarry Smith { 70421e72a00SBarry Smith PetscErrorCode ierr; 70521e72a00SBarry Smith PetscInt *L_ = *L, ak, bk, k; 70621e72a00SBarry Smith 70721e72a00SBarry Smith if (!L_) { 70821e72a00SBarry Smith ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr); 70921e72a00SBarry Smith L_ = *L; 71021e72a00SBarry Smith } 71121e72a00SBarry Smith k = ak = bk = 0; 71221e72a00SBarry Smith while (ak < an && bk < bn) { 71321e72a00SBarry Smith if (aI[ak] == bI[bk]) { 71421e72a00SBarry Smith L_[k] = aI[ak]; 71521e72a00SBarry Smith ++ak; 71621e72a00SBarry Smith ++bk; 71721e72a00SBarry Smith ++k; 71821e72a00SBarry Smith } else if (aI[ak] < bI[bk]) { 71921e72a00SBarry Smith L_[k] = aI[ak]; 72021e72a00SBarry Smith ++ak; 72121e72a00SBarry Smith ++k; 72221e72a00SBarry Smith } else { 72321e72a00SBarry Smith L_[k] = bI[bk]; 72421e72a00SBarry Smith ++bk; 72521e72a00SBarry Smith ++k; 72621e72a00SBarry Smith } 72721e72a00SBarry Smith } 72821e72a00SBarry Smith if (ak < an) { 72921e72a00SBarry Smith ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr); 73021e72a00SBarry Smith k += (an-ak); 73121e72a00SBarry Smith } 73221e72a00SBarry Smith if (bk < bn) { 73321e72a00SBarry Smith ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr); 73421e72a00SBarry Smith k += (bn-bk); 73521e72a00SBarry Smith } 73621e72a00SBarry Smith *n = k; 73721e72a00SBarry Smith PetscFunctionReturn(0); 73821e72a00SBarry Smith } 739b4301105SBarry Smith 740c1f0200aSDmitry Karpeev /*@ 74121e72a00SBarry Smith PetscMergeIntArrayPair - Merges two SORTED integer arrays that share NO common values along with an additional array of integers. 742c1f0200aSDmitry Karpeev The additional arrays are the same length as sorted arrays and are merged 743c1f0200aSDmitry Karpeev in the order determined by the merging of the sorted pair. 744c1f0200aSDmitry Karpeev 745c1f0200aSDmitry Karpeev Not Collective 746c1f0200aSDmitry Karpeev 747c1f0200aSDmitry Karpeev Input Parameters: 748c1f0200aSDmitry Karpeev + an - number of values in the first array 749c1f0200aSDmitry Karpeev . aI - first sorted array of integers 750c1f0200aSDmitry Karpeev . aJ - first additional array of integers 751c1f0200aSDmitry Karpeev . bn - number of values in the second array 752c1f0200aSDmitry Karpeev . bI - second array of integers 753c1f0200aSDmitry Karpeev - bJ - second additional of integers 754c1f0200aSDmitry Karpeev 755c1f0200aSDmitry Karpeev Output Parameters: 756c1f0200aSDmitry Karpeev + n - number of values in the merged array (== an + bn) 75714c49006SJed Brown . L - merged sorted array 758c1f0200aSDmitry Karpeev - J - merged additional array 759c1f0200aSDmitry Karpeev 76070d8d27cSBarry Smith 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 761c1f0200aSDmitry Karpeev Level: intermediate 762c1f0200aSDmitry Karpeev 763c1f0200aSDmitry Karpeev Concepts: merging^arrays 764c1f0200aSDmitry Karpeev 765c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 766c1f0200aSDmitry Karpeev @*/ 7676c2863d0SBarry 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) 768c1f0200aSDmitry Karpeev { 769c1f0200aSDmitry Karpeev PetscErrorCode ierr; 77070d8d27cSBarry Smith PetscInt n_, *L_, *J_, ak, bk, k; 771c1f0200aSDmitry Karpeev 77214c49006SJed Brown PetscFunctionBegin; 77370d8d27cSBarry Smith PetscValidIntPointer(L,8); 77470d8d27cSBarry Smith PetscValidIntPointer(J,9); 775c1f0200aSDmitry Karpeev n_ = an + bn; 776c1f0200aSDmitry Karpeev *n = n_; 77770d8d27cSBarry Smith if (!*L) { 778785e854fSJed Brown ierr = PetscMalloc1(n_, L);CHKERRQ(ierr); 77970d8d27cSBarry Smith } 780d7aa01a8SHong Zhang L_ = *L; 78170d8d27cSBarry Smith if (!*J) { 78270d8d27cSBarry Smith ierr = PetscMalloc1(n_, J);CHKERRQ(ierr); 783c1f0200aSDmitry Karpeev } 784c1f0200aSDmitry Karpeev J_ = *J; 785c1f0200aSDmitry Karpeev k = ak = bk = 0; 786c1f0200aSDmitry Karpeev while (ak < an && bk < bn) { 787c1f0200aSDmitry Karpeev if (aI[ak] <= bI[bk]) { 788d7aa01a8SHong Zhang L_[k] = aI[ak]; 789c1f0200aSDmitry Karpeev J_[k] = aJ[ak]; 790c1f0200aSDmitry Karpeev ++ak; 791c1f0200aSDmitry Karpeev ++k; 792a297a907SKarl Rupp } else { 793d7aa01a8SHong Zhang L_[k] = bI[bk]; 794c1f0200aSDmitry Karpeev J_[k] = bJ[bk]; 795c1f0200aSDmitry Karpeev ++bk; 796c1f0200aSDmitry Karpeev ++k; 797c1f0200aSDmitry Karpeev } 798c1f0200aSDmitry Karpeev } 799c1f0200aSDmitry Karpeev if (ak < an) { 800d7aa01a8SHong Zhang ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr); 801c1f0200aSDmitry Karpeev ierr = PetscMemcpy(J_+k,aJ+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr); 802c1f0200aSDmitry Karpeev k += (an-ak); 803c1f0200aSDmitry Karpeev } 804c1f0200aSDmitry Karpeev if (bk < bn) { 805d7aa01a8SHong Zhang ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr); 806c1f0200aSDmitry Karpeev ierr = PetscMemcpy(J_+k,bJ+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr); 807c1f0200aSDmitry Karpeev } 808c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 809c1f0200aSDmitry Karpeev } 810c1f0200aSDmitry Karpeev 811e498c390SJed Brown /*@ 812e498c390SJed Brown PetscMergeMPIIntArray - Merges two SORTED integer arrays. 813e498c390SJed Brown 814e498c390SJed Brown Not Collective 815e498c390SJed Brown 816e498c390SJed Brown Input Parameters: 817e498c390SJed Brown + an - number of values in the first array 818e498c390SJed Brown . aI - first sorted array of integers 819e498c390SJed Brown . bn - number of values in the second array 820e498c390SJed Brown - bI - second array of integers 821e498c390SJed Brown 822e498c390SJed Brown Output Parameters: 823e498c390SJed Brown + n - number of values in the merged array (<= an + bn) 824e498c390SJed Brown - L - merged sorted array, allocated if address of NULL pointer is passed 825e498c390SJed Brown 826e498c390SJed Brown Level: intermediate 827e498c390SJed Brown 828e498c390SJed Brown Concepts: merging^arrays 829e498c390SJed Brown 830e498c390SJed Brown .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 831e498c390SJed Brown @*/ 832e498c390SJed Brown PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L) 833e498c390SJed Brown { 834e498c390SJed Brown PetscErrorCode ierr; 835e498c390SJed Brown PetscInt ai,bi,k; 836e498c390SJed Brown 837e498c390SJed Brown PetscFunctionBegin; 838e498c390SJed Brown if (!*L) {ierr = PetscMalloc1((an+bn),L);CHKERRQ(ierr);} 839e498c390SJed Brown for (ai=0,bi=0,k=0; ai<an || bi<bn; ) { 840e498c390SJed Brown PetscInt t = -1; 841e498c390SJed Brown for ( ; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai]; 842e498c390SJed Brown for ( ; bi<bn && bI[bi] == t; bi++); 843e498c390SJed Brown for ( ; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi]; 844e498c390SJed Brown for ( ; ai<an && aI[ai] == t; ai++); 845e498c390SJed Brown } 846e498c390SJed Brown *n = k; 847e498c390SJed Brown PetscFunctionReturn(0); 848e498c390SJed Brown } 849e5c89e4eSSatish Balay 8506c2863d0SBarry Smith /*@C 85142eaa7f3SBarry Smith PetscProcessTree - Prepares tree data to be displayed graphically 85242eaa7f3SBarry Smith 85342eaa7f3SBarry Smith Not Collective 85442eaa7f3SBarry Smith 85542eaa7f3SBarry Smith Input Parameters: 85642eaa7f3SBarry Smith + n - number of values 85742eaa7f3SBarry Smith . mask - indicates those entries in the tree, location 0 is always masked 85842eaa7f3SBarry Smith - parentid - indicates the parent of each entry 85942eaa7f3SBarry Smith 86042eaa7f3SBarry Smith Output Parameters: 86142eaa7f3SBarry Smith + Nlevels - the number of levels 86242eaa7f3SBarry Smith . Level - for each node tells its level 86342eaa7f3SBarry Smith . Levelcnts - the number of nodes on each level 86442eaa7f3SBarry Smith . Idbylevel - a list of ids on each of the levels, first level followed by second etc 86542eaa7f3SBarry Smith - Column - for each id tells its column index 86642eaa7f3SBarry Smith 8676c2863d0SBarry Smith Level: developer 86842eaa7f3SBarry Smith 8696c2863d0SBarry Smith Notes: This code is not currently used 87042eaa7f3SBarry Smith 87142eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation() 87242eaa7f3SBarry Smith @*/ 8737087cfbeSBarry Smith PetscErrorCode PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column) 87442eaa7f3SBarry Smith { 87542eaa7f3SBarry Smith PetscInt i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column; 87642eaa7f3SBarry Smith PetscErrorCode ierr; 877ace3abfcSBarry Smith PetscBool done = PETSC_FALSE; 87842eaa7f3SBarry Smith 87942eaa7f3SBarry Smith PetscFunctionBegin; 88042eaa7f3SBarry Smith if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set"); 88142eaa7f3SBarry Smith for (i=0; i<n; i++) { 88242eaa7f3SBarry Smith if (mask[i]) continue; 88342eaa7f3SBarry Smith if (parentid[i] == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent"); 88442eaa7f3SBarry Smith if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked"); 88542eaa7f3SBarry Smith } 88642eaa7f3SBarry Smith 88742eaa7f3SBarry Smith for (i=0; i<n; i++) { 88842eaa7f3SBarry Smith if (!mask[i]) nmask++; 88942eaa7f3SBarry Smith } 89042eaa7f3SBarry Smith 89142eaa7f3SBarry Smith /* determine the level in the tree of each node */ 8921795a4d1SJed Brown ierr = PetscCalloc1(n,&level);CHKERRQ(ierr); 893a297a907SKarl Rupp 89442eaa7f3SBarry Smith level[0] = 1; 89542eaa7f3SBarry Smith while (!done) { 89642eaa7f3SBarry Smith done = PETSC_TRUE; 89742eaa7f3SBarry Smith for (i=0; i<n; i++) { 89842eaa7f3SBarry Smith if (mask[i]) continue; 89942eaa7f3SBarry Smith if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1; 90042eaa7f3SBarry Smith else if (!level[i]) done = PETSC_FALSE; 90142eaa7f3SBarry Smith } 90242eaa7f3SBarry Smith } 90342eaa7f3SBarry Smith for (i=0; i<n; i++) { 90442eaa7f3SBarry Smith level[i]--; 90542eaa7f3SBarry Smith nlevels = PetscMax(nlevels,level[i]); 90642eaa7f3SBarry Smith } 90742eaa7f3SBarry Smith 90842eaa7f3SBarry Smith /* count the number of nodes on each level and its max */ 9091795a4d1SJed Brown ierr = PetscCalloc1(nlevels,&levelcnt);CHKERRQ(ierr); 91042eaa7f3SBarry Smith for (i=0; i<n; i++) { 91142eaa7f3SBarry Smith if (mask[i]) continue; 91242eaa7f3SBarry Smith levelcnt[level[i]-1]++; 91342eaa7f3SBarry Smith } 914a297a907SKarl Rupp for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]); 91542eaa7f3SBarry Smith 91642eaa7f3SBarry Smith /* for each level sort the ids by the parent id */ 917dcca6d9dSJed Brown ierr = PetscMalloc2(levelmax,&workid,levelmax,&workparentid);CHKERRQ(ierr); 918785e854fSJed Brown ierr = PetscMalloc1(nmask,&idbylevel);CHKERRQ(ierr); 91942eaa7f3SBarry Smith for (j=1; j<=nlevels;j++) { 92042eaa7f3SBarry Smith cnt = 0; 92142eaa7f3SBarry Smith for (i=0; i<n; i++) { 92242eaa7f3SBarry Smith if (mask[i]) continue; 92342eaa7f3SBarry Smith if (level[i] != j) continue; 92442eaa7f3SBarry Smith workid[cnt] = i; 92542eaa7f3SBarry Smith workparentid[cnt++] = parentid[i]; 92642eaa7f3SBarry Smith } 92742eaa7f3SBarry Smith /* PetscIntView(cnt,workparentid,0); 92842eaa7f3SBarry Smith PetscIntView(cnt,workid,0); 92942eaa7f3SBarry Smith ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr); 93042eaa7f3SBarry Smith PetscIntView(cnt,workparentid,0); 93142eaa7f3SBarry Smith PetscIntView(cnt,workid,0);*/ 93242eaa7f3SBarry Smith ierr = PetscMemcpy(idbylevel+tcnt,workid,cnt*sizeof(PetscInt));CHKERRQ(ierr); 93342eaa7f3SBarry Smith tcnt += cnt; 93442eaa7f3SBarry Smith } 93542eaa7f3SBarry Smith if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes"); 93642eaa7f3SBarry Smith ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr); 93742eaa7f3SBarry Smith 93842eaa7f3SBarry Smith /* for each node list its column */ 939785e854fSJed Brown ierr = PetscMalloc1(n,&column);CHKERRQ(ierr); 94042eaa7f3SBarry Smith cnt = 0; 94142eaa7f3SBarry Smith for (j=0; j<nlevels; j++) { 94242eaa7f3SBarry Smith for (i=0; i<levelcnt[j]; i++) { 94342eaa7f3SBarry Smith column[idbylevel[cnt++]] = i; 94442eaa7f3SBarry Smith } 94542eaa7f3SBarry Smith } 94642eaa7f3SBarry Smith 94742eaa7f3SBarry Smith *Nlevels = nlevels; 94842eaa7f3SBarry Smith *Level = level; 94942eaa7f3SBarry Smith *Levelcnt = levelcnt; 95042eaa7f3SBarry Smith *Idbylevel = idbylevel; 95142eaa7f3SBarry Smith *Column = column; 95242eaa7f3SBarry Smith PetscFunctionReturn(0); 95342eaa7f3SBarry Smith } 954