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 /*@ 861db36a52SBarry Smith PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries 871db36a52SBarry Smith 881db36a52SBarry Smith Not Collective 891db36a52SBarry Smith 901db36a52SBarry Smith Input Parameters: 911db36a52SBarry Smith + n - number of values 921db36a52SBarry Smith - ii - array of integers 931db36a52SBarry Smith 941db36a52SBarry Smith Output Parameter: 951db36a52SBarry Smith . n - number of non-redundant values 961db36a52SBarry Smith 971db36a52SBarry Smith Level: intermediate 981db36a52SBarry Smith 991db36a52SBarry Smith Concepts: sorting^ints 1001db36a52SBarry Smith 1011db36a52SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt() 1021db36a52SBarry Smith @*/ 1037087cfbeSBarry Smith PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n,PetscInt ii[]) 1041db36a52SBarry Smith { 1051db36a52SBarry Smith PetscErrorCode ierr; 1061db36a52SBarry Smith PetscInt i,s = 0,N = *n, b = 0; 1071db36a52SBarry Smith 1081db36a52SBarry Smith PetscFunctionBegin; 1091db36a52SBarry Smith ierr = PetscSortInt(N,ii);CHKERRQ(ierr); 1101db36a52SBarry Smith for (i=0; i<N-1; i++) { 111a297a907SKarl Rupp if (ii[b+s+1] != ii[b]) { 112a297a907SKarl Rupp ii[b+1] = ii[b+s+1]; b++; 113a297a907SKarl Rupp } else s++; 1141db36a52SBarry Smith } 1151db36a52SBarry Smith *n = N - s; 1161db36a52SBarry Smith PetscFunctionReturn(0); 1171db36a52SBarry Smith } 1181db36a52SBarry Smith 11960e03357SMatthew G Knepley /*@ 120d9f1162dSJed Brown PetscFindInt - Finds integer in a sorted array of integers 12160e03357SMatthew G Knepley 12260e03357SMatthew G Knepley Not Collective 12360e03357SMatthew G Knepley 12460e03357SMatthew G Knepley Input Parameters: 12560e03357SMatthew G Knepley + key - the integer to locate 126d9f1162dSJed Brown . n - number of values in the array 12760e03357SMatthew G Knepley - ii - array of integers 12860e03357SMatthew G Knepley 12960e03357SMatthew G Knepley Output Parameter: 130d9f1162dSJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 13160e03357SMatthew G Knepley 13260e03357SMatthew G Knepley Level: intermediate 13360e03357SMatthew G Knepley 13460e03357SMatthew G Knepley Concepts: sorting^ints 13560e03357SMatthew G Knepley 13660e03357SMatthew G Knepley .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt() 13760e03357SMatthew G Knepley @*/ 138026ec6cbSMatthew G Knepley PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt ii[], PetscInt *loc) 13960e03357SMatthew G Knepley { 140d9f1162dSJed Brown PetscInt lo = 0,hi = n; 14160e03357SMatthew G Knepley 14260e03357SMatthew G Knepley PetscFunctionBegin; 14360e03357SMatthew G Knepley PetscValidPointer(loc,4); 14498781d82SMatthew G Knepley if (!n) {*loc = -1; PetscFunctionReturn(0);} 14598781d82SMatthew G Knepley PetscValidPointer(ii,3); 146d9f1162dSJed Brown while (hi - lo > 1) { 147d9f1162dSJed Brown PetscInt mid = lo + (hi - lo)/2; 148d9f1162dSJed Brown if (key < ii[mid]) hi = mid; 149d9f1162dSJed Brown else lo = mid; 15060e03357SMatthew G Knepley } 151d9f1162dSJed Brown *loc = key == ii[lo] ? lo : -(lo + (key > ii[lo]) + 1); 15260e03357SMatthew G Knepley PetscFunctionReturn(0); 15360e03357SMatthew G Knepley } 15460e03357SMatthew G Knepley 155*d2aeb606SJed Brown #undef __FUNCT__ 156*d2aeb606SJed Brown #define __FUNCT__ "PetscFindMPIInt" 157*d2aeb606SJed Brown /*@ 158*d2aeb606SJed Brown PetscFindMPIInt - Finds MPI integer in a sorted array of integers 159*d2aeb606SJed Brown 160*d2aeb606SJed Brown Not Collective 161*d2aeb606SJed Brown 162*d2aeb606SJed Brown Input Parameters: 163*d2aeb606SJed Brown + key - the integer to locate 164*d2aeb606SJed Brown . n - number of values in the array 165*d2aeb606SJed Brown - ii - array of integers 166*d2aeb606SJed Brown 167*d2aeb606SJed Brown Output Parameter: 168*d2aeb606SJed Brown . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 169*d2aeb606SJed Brown 170*d2aeb606SJed Brown Level: intermediate 171*d2aeb606SJed Brown 172*d2aeb606SJed Brown Concepts: sorting^ints 173*d2aeb606SJed Brown 174*d2aeb606SJed Brown .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt() 175*d2aeb606SJed Brown @*/ 176*d2aeb606SJed Brown PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt ii[], PetscInt *loc) 177*d2aeb606SJed Brown { 178*d2aeb606SJed Brown PetscInt lo = 0,hi = n; 179*d2aeb606SJed Brown 180*d2aeb606SJed Brown PetscFunctionBegin; 181*d2aeb606SJed Brown PetscValidPointer(loc,4); 182*d2aeb606SJed Brown if (!n) {*loc = -1; PetscFunctionReturn(0);} 183*d2aeb606SJed Brown PetscValidPointer(ii,3); 184*d2aeb606SJed Brown while (hi - lo > 1) { 185*d2aeb606SJed Brown PetscInt mid = lo + (hi - lo)/2; 186*d2aeb606SJed Brown if (key < ii[mid]) hi = mid; 187*d2aeb606SJed Brown else lo = mid; 188*d2aeb606SJed Brown } 189*d2aeb606SJed Brown *loc = key == ii[lo] ? lo : -(lo + (key > ii[lo]) + 1); 190*d2aeb606SJed Brown PetscFunctionReturn(0); 191*d2aeb606SJed Brown } 1921db36a52SBarry Smith 193e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 194e5c89e4eSSatish Balay #define SWAP2(a,b,c,d,t) {t=a;a=b;b=t;t=c;c=d;d=t;} 195e5c89e4eSSatish Balay 196e5c89e4eSSatish Balay /* 197e5c89e4eSSatish Balay A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 198e5c89e4eSSatish Balay Assumes 0 origin for v, number of elements = right+1 (right is index of 199e5c89e4eSSatish Balay right-most member). 200e5c89e4eSSatish Balay */ 201e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithArray_Private(PetscInt *v,PetscInt *V,PetscInt right) 202e5c89e4eSSatish Balay { 203e5c89e4eSSatish Balay PetscErrorCode ierr; 204e5c89e4eSSatish Balay PetscInt i,vl,last,tmp; 205e5c89e4eSSatish Balay 206e5c89e4eSSatish Balay PetscFunctionBegin; 207e5c89e4eSSatish Balay if (right <= 1) { 208e5c89e4eSSatish Balay if (right == 1) { 209e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp); 210e5c89e4eSSatish Balay } 211e5c89e4eSSatish Balay PetscFunctionReturn(0); 212e5c89e4eSSatish Balay } 213e5c89e4eSSatish Balay SWAP2(v[0],v[right/2],V[0],V[right/2],tmp); 214e5c89e4eSSatish Balay vl = v[0]; 215e5c89e4eSSatish Balay last = 0; 216e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 217e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);} 218e5c89e4eSSatish Balay } 219e5c89e4eSSatish Balay SWAP2(v[0],v[last],V[0],V[last],tmp); 220e5c89e4eSSatish Balay ierr = PetscSortIntWithArray_Private(v,V,last-1);CHKERRQ(ierr); 221e5c89e4eSSatish Balay ierr = PetscSortIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 222e5c89e4eSSatish Balay PetscFunctionReturn(0); 223e5c89e4eSSatish Balay } 224e5c89e4eSSatish Balay 225e5c89e4eSSatish Balay /*@ 226e5c89e4eSSatish Balay PetscSortIntWithArray - Sorts an array of integers in place in increasing order; 227e5c89e4eSSatish Balay changes a second array to match the sorted first array. 228e5c89e4eSSatish Balay 229e5c89e4eSSatish Balay Not Collective 230e5c89e4eSSatish Balay 231e5c89e4eSSatish Balay Input Parameters: 232e5c89e4eSSatish Balay + n - number of values 233e5c89e4eSSatish Balay . i - array of integers 234e5c89e4eSSatish Balay - I - second array of integers 235e5c89e4eSSatish Balay 236e5c89e4eSSatish Balay Level: intermediate 237e5c89e4eSSatish Balay 238e5c89e4eSSatish Balay Concepts: sorting^ints with array 239e5c89e4eSSatish Balay 240e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 241e5c89e4eSSatish Balay @*/ 2427087cfbeSBarry Smith PetscErrorCode PetscSortIntWithArray(PetscInt n,PetscInt i[],PetscInt Ii[]) 243e5c89e4eSSatish Balay { 244e5c89e4eSSatish Balay PetscErrorCode ierr; 245e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 246e5c89e4eSSatish Balay 247e5c89e4eSSatish Balay PetscFunctionBegin; 248e5c89e4eSSatish Balay if (n<8) { 249e5c89e4eSSatish Balay for (k=0; k<n; k++) { 250e5c89e4eSSatish Balay ik = i[k]; 251e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 252e5c89e4eSSatish Balay if (ik > i[j]) { 253b7940d39SSatish Balay SWAP2(i[k],i[j],Ii[k],Ii[j],tmp); 254e5c89e4eSSatish Balay ik = i[k]; 255e5c89e4eSSatish Balay } 256e5c89e4eSSatish Balay } 257e5c89e4eSSatish Balay } 258e5c89e4eSSatish Balay } else { 259b7940d39SSatish Balay ierr = PetscSortIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr); 260e5c89e4eSSatish Balay } 261e5c89e4eSSatish Balay PetscFunctionReturn(0); 262e5c89e4eSSatish Balay } 263e5c89e4eSSatish Balay 264c1f0200aSDmitry Karpeev 265c1f0200aSDmitry 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;} 266c1f0200aSDmitry Karpeev 267c1f0200aSDmitry Karpeev /* 268c1f0200aSDmitry Karpeev A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 269c1f0200aSDmitry Karpeev Assumes 0 origin for v, number of elements = right+1 (right is index of 270c1f0200aSDmitry Karpeev right-most member). 271c1f0200aSDmitry Karpeev */ 272d7aa01a8SHong Zhang static PetscErrorCode PetscSortIntWithArrayPair_Private(PetscInt *L,PetscInt *J, PetscInt *K, PetscInt right) 273c1f0200aSDmitry Karpeev { 274c1f0200aSDmitry Karpeev PetscErrorCode ierr; 275c1f0200aSDmitry Karpeev PetscInt i,vl,last,tmp; 276c1f0200aSDmitry Karpeev 277c1f0200aSDmitry Karpeev PetscFunctionBegin; 278c1f0200aSDmitry Karpeev if (right <= 1) { 279c1f0200aSDmitry Karpeev if (right == 1) { 280d7aa01a8SHong Zhang if (L[0] > L[1]) SWAP3(L[0],L[1],J[0],J[1],K[0],K[1],tmp); 281c1f0200aSDmitry Karpeev } 282c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 283c1f0200aSDmitry Karpeev } 284d7aa01a8SHong Zhang SWAP3(L[0],L[right/2],J[0],J[right/2],K[0],K[right/2],tmp); 285d7aa01a8SHong Zhang vl = L[0]; 286c1f0200aSDmitry Karpeev last = 0; 287c1f0200aSDmitry Karpeev for (i=1; i<=right; i++) { 288d7aa01a8SHong Zhang if (L[i] < vl) {last++; SWAP3(L[last],L[i],J[last],J[i],K[last],K[i],tmp);} 289c1f0200aSDmitry Karpeev } 290d7aa01a8SHong Zhang SWAP3(L[0],L[last],J[0],J[last],K[0],K[last],tmp); 291d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L,J,K,last-1);CHKERRQ(ierr); 292d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L+last+1,J+last+1,K+last+1,right-(last+1));CHKERRQ(ierr); 293c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 294c1f0200aSDmitry Karpeev } 295c1f0200aSDmitry Karpeev 296c1f0200aSDmitry Karpeev /*@ 297c1f0200aSDmitry Karpeev PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order; 298c1f0200aSDmitry Karpeev changes a pair of integer arrays to match the sorted first array. 299c1f0200aSDmitry Karpeev 300c1f0200aSDmitry Karpeev Not Collective 301c1f0200aSDmitry Karpeev 302c1f0200aSDmitry Karpeev Input Parameters: 303c1f0200aSDmitry Karpeev + n - number of values 304c1f0200aSDmitry Karpeev . I - array of integers 305c1f0200aSDmitry Karpeev . J - second array of integers (first array of the pair) 306c1f0200aSDmitry Karpeev - K - third array of integers (second array of the pair) 307c1f0200aSDmitry Karpeev 308c1f0200aSDmitry Karpeev Level: intermediate 309c1f0200aSDmitry Karpeev 310c1f0200aSDmitry Karpeev Concepts: sorting^ints with array pair 311c1f0200aSDmitry Karpeev 312c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray() 313c1f0200aSDmitry Karpeev @*/ 3146c2863d0SBarry Smith PetscErrorCode PetscSortIntWithArrayPair(PetscInt n,PetscInt L[],PetscInt J[], PetscInt K[]) 315c1f0200aSDmitry Karpeev { 316c1f0200aSDmitry Karpeev PetscErrorCode ierr; 317c1f0200aSDmitry Karpeev PetscInt j,k,tmp,ik; 318c1f0200aSDmitry Karpeev 319c1f0200aSDmitry Karpeev PetscFunctionBegin; 320c1f0200aSDmitry Karpeev if (n<8) { 321c1f0200aSDmitry Karpeev for (k=0; k<n; k++) { 322d7aa01a8SHong Zhang ik = L[k]; 323c1f0200aSDmitry Karpeev for (j=k+1; j<n; j++) { 324d7aa01a8SHong Zhang if (ik > L[j]) { 325d7aa01a8SHong Zhang SWAP3(L[k],L[j],J[k],J[j],K[k],K[j],tmp); 326d7aa01a8SHong Zhang ik = L[k]; 327c1f0200aSDmitry Karpeev } 328c1f0200aSDmitry Karpeev } 329c1f0200aSDmitry Karpeev } 330c1f0200aSDmitry Karpeev } else { 331d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L,J,K,n-1);CHKERRQ(ierr); 332c1f0200aSDmitry Karpeev } 333c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 334c1f0200aSDmitry Karpeev } 335c1f0200aSDmitry Karpeev 33617d7d925SStefano Zampini /* 33717d7d925SStefano Zampini A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 33817d7d925SStefano Zampini Assumes 0 origin for v, number of elements = right+1 (right is index of 33917d7d925SStefano Zampini right-most member). 34017d7d925SStefano Zampini */ 34117d7d925SStefano Zampini static void PetscSortMPIInt_Private(PetscMPIInt *v,PetscInt right) 34217d7d925SStefano Zampini { 34317d7d925SStefano Zampini PetscInt i,j; 34417d7d925SStefano Zampini PetscMPIInt pivot,tmp; 34517d7d925SStefano Zampini 34617d7d925SStefano Zampini if (right <= 1) { 34717d7d925SStefano Zampini if (right == 1) { 34817d7d925SStefano Zampini if (v[0] > v[1]) SWAP(v[0],v[1],tmp); 34917d7d925SStefano Zampini } 35017d7d925SStefano Zampini return; 35117d7d925SStefano Zampini } 35217d7d925SStefano Zampini i = MEDIAN(v,right); /* Choose a pivot */ 35317d7d925SStefano Zampini SWAP(v[0],v[i],tmp); /* Move it out of the way */ 35417d7d925SStefano Zampini pivot = v[0]; 35517d7d925SStefano Zampini for (i=0,j=right+1;; ) { 35617d7d925SStefano Zampini while (++i < j && v[i] <= pivot) ; /* Scan from the left */ 35717d7d925SStefano Zampini while (v[--j] > pivot) ; /* Scan from the right */ 35817d7d925SStefano Zampini if (i >= j) break; 35917d7d925SStefano Zampini SWAP(v[i],v[j],tmp); 36017d7d925SStefano Zampini } 36117d7d925SStefano Zampini SWAP(v[0],v[j],tmp); /* Put pivot back in place. */ 36217d7d925SStefano Zampini PetscSortMPIInt_Private(v,j-1); 36317d7d925SStefano Zampini PetscSortMPIInt_Private(v+j+1,right-(j+1)); 36417d7d925SStefano Zampini } 36517d7d925SStefano Zampini 36617d7d925SStefano Zampini /*@ 36717d7d925SStefano Zampini PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order. 36817d7d925SStefano Zampini 36917d7d925SStefano Zampini Not Collective 37017d7d925SStefano Zampini 37117d7d925SStefano Zampini Input Parameters: 37217d7d925SStefano Zampini + n - number of values 37317d7d925SStefano Zampini - i - array of integers 37417d7d925SStefano Zampini 37517d7d925SStefano Zampini Level: intermediate 37617d7d925SStefano Zampini 37717d7d925SStefano Zampini Concepts: sorting^ints 37817d7d925SStefano Zampini 37917d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation() 38017d7d925SStefano Zampini @*/ 38117d7d925SStefano Zampini PetscErrorCode PetscSortMPIInt(PetscInt n,PetscMPIInt i[]) 38217d7d925SStefano Zampini { 38317d7d925SStefano Zampini PetscInt j,k; 38417d7d925SStefano Zampini PetscMPIInt tmp,ik; 38517d7d925SStefano Zampini 38617d7d925SStefano Zampini PetscFunctionBegin; 38717d7d925SStefano Zampini if (n<8) { 38817d7d925SStefano Zampini for (k=0; k<n; k++) { 38917d7d925SStefano Zampini ik = i[k]; 39017d7d925SStefano Zampini for (j=k+1; j<n; j++) { 39117d7d925SStefano Zampini if (ik > i[j]) { 39217d7d925SStefano Zampini SWAP(i[k],i[j],tmp); 39317d7d925SStefano Zampini ik = i[k]; 39417d7d925SStefano Zampini } 39517d7d925SStefano Zampini } 39617d7d925SStefano Zampini } 397a297a907SKarl Rupp } else PetscSortMPIInt_Private(i,n-1); 39817d7d925SStefano Zampini PetscFunctionReturn(0); 39917d7d925SStefano Zampini } 40017d7d925SStefano Zampini 40117d7d925SStefano Zampini /*@ 40217d7d925SStefano Zampini PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries 40317d7d925SStefano Zampini 40417d7d925SStefano Zampini Not Collective 40517d7d925SStefano Zampini 40617d7d925SStefano Zampini Input Parameters: 40717d7d925SStefano Zampini + n - number of values 40817d7d925SStefano Zampini - ii - array of integers 40917d7d925SStefano Zampini 41017d7d925SStefano Zampini Output Parameter: 41117d7d925SStefano Zampini . n - number of non-redundant values 41217d7d925SStefano Zampini 41317d7d925SStefano Zampini Level: intermediate 41417d7d925SStefano Zampini 41517d7d925SStefano Zampini Concepts: sorting^ints 41617d7d925SStefano Zampini 41717d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt() 41817d7d925SStefano Zampini @*/ 41917d7d925SStefano Zampini PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt ii[]) 42017d7d925SStefano Zampini { 42117d7d925SStefano Zampini PetscErrorCode ierr; 42217d7d925SStefano Zampini PetscInt i,s = 0,N = *n, b = 0; 42317d7d925SStefano Zampini 42417d7d925SStefano Zampini PetscFunctionBegin; 42517d7d925SStefano Zampini ierr = PetscSortMPIInt(N,ii);CHKERRQ(ierr); 42617d7d925SStefano Zampini for (i=0; i<N-1; i++) { 427a297a907SKarl Rupp if (ii[b+s+1] != ii[b]) { 428a297a907SKarl Rupp ii[b+1] = ii[b+s+1]; b++; 429a297a907SKarl Rupp } else s++; 43017d7d925SStefano Zampini } 43117d7d925SStefano Zampini *n = N - s; 43217d7d925SStefano Zampini PetscFunctionReturn(0); 43317d7d925SStefano Zampini } 434c1f0200aSDmitry Karpeev 4354d615ea0SBarry Smith /* 4364d615ea0SBarry Smith A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 4374d615ea0SBarry Smith Assumes 0 origin for v, number of elements = right+1 (right is index of 4384d615ea0SBarry Smith right-most member). 4394d615ea0SBarry Smith */ 4404d615ea0SBarry Smith static PetscErrorCode PetscSortMPIIntWithArray_Private(PetscMPIInt *v,PetscMPIInt *V,PetscMPIInt right) 4414d615ea0SBarry Smith { 4424d615ea0SBarry Smith PetscErrorCode ierr; 4434d615ea0SBarry Smith PetscMPIInt i,vl,last,tmp; 4444d615ea0SBarry Smith 4454d615ea0SBarry Smith PetscFunctionBegin; 4464d615ea0SBarry Smith if (right <= 1) { 4474d615ea0SBarry Smith if (right == 1) { 4484d615ea0SBarry Smith if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp); 4494d615ea0SBarry Smith } 4504d615ea0SBarry Smith PetscFunctionReturn(0); 4514d615ea0SBarry Smith } 4524d615ea0SBarry Smith SWAP2(v[0],v[right/2],V[0],V[right/2],tmp); 4534d615ea0SBarry Smith vl = v[0]; 4544d615ea0SBarry Smith last = 0; 4554d615ea0SBarry Smith for (i=1; i<=right; i++) { 4564d615ea0SBarry Smith if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);} 4574d615ea0SBarry Smith } 4584d615ea0SBarry Smith SWAP2(v[0],v[last],V[0],V[last],tmp); 4594d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v,V,last-1);CHKERRQ(ierr); 4604d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 4614d615ea0SBarry Smith PetscFunctionReturn(0); 4624d615ea0SBarry Smith } 4634d615ea0SBarry Smith 4644d615ea0SBarry Smith /*@ 4654d615ea0SBarry Smith PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order; 4664d615ea0SBarry Smith changes a second array to match the sorted first array. 4674d615ea0SBarry Smith 4684d615ea0SBarry Smith Not Collective 4694d615ea0SBarry Smith 4704d615ea0SBarry Smith Input Parameters: 4714d615ea0SBarry Smith + n - number of values 4724d615ea0SBarry Smith . i - array of integers 4734d615ea0SBarry Smith - I - second array of integers 4744d615ea0SBarry Smith 4754d615ea0SBarry Smith Level: intermediate 4764d615ea0SBarry Smith 4774d615ea0SBarry Smith Concepts: sorting^ints with array 4784d615ea0SBarry Smith 4794d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 4804d615ea0SBarry Smith @*/ 4817087cfbeSBarry Smith PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt i[],PetscMPIInt Ii[]) 4824d615ea0SBarry Smith { 4834d615ea0SBarry Smith PetscErrorCode ierr; 4844d615ea0SBarry Smith PetscMPIInt j,k,tmp,ik; 4854d615ea0SBarry Smith 4864d615ea0SBarry Smith PetscFunctionBegin; 4874d615ea0SBarry Smith if (n<8) { 4884d615ea0SBarry Smith for (k=0; k<n; k++) { 4894d615ea0SBarry Smith ik = i[k]; 4904d615ea0SBarry Smith for (j=k+1; j<n; j++) { 4914d615ea0SBarry Smith if (ik > i[j]) { 4924d615ea0SBarry Smith SWAP2(i[k],i[j],Ii[k],Ii[j],tmp); 4934d615ea0SBarry Smith ik = i[k]; 4944d615ea0SBarry Smith } 4954d615ea0SBarry Smith } 4964d615ea0SBarry Smith } 4974d615ea0SBarry Smith } else { 4984d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr); 4994d615ea0SBarry Smith } 5004d615ea0SBarry Smith PetscFunctionReturn(0); 5014d615ea0SBarry Smith } 5024d615ea0SBarry Smith 503e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 504e5c89e4eSSatish Balay #define SWAP2IntScalar(a,b,c,d,t,ts) {t=a;a=b;b=t;ts=c;c=d;d=ts;} 505e5c89e4eSSatish Balay 506e5c89e4eSSatish Balay /* 507e5c89e4eSSatish Balay Modified from PetscSortIntWithArray_Private(). 508e5c89e4eSSatish Balay */ 509e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithScalarArray_Private(PetscInt *v,PetscScalar *V,PetscInt right) 510e5c89e4eSSatish Balay { 511e5c89e4eSSatish Balay PetscErrorCode ierr; 512e5c89e4eSSatish Balay PetscInt i,vl,last,tmp; 513e5c89e4eSSatish Balay PetscScalar stmp; 514e5c89e4eSSatish Balay 515e5c89e4eSSatish Balay PetscFunctionBegin; 516e5c89e4eSSatish Balay if (right <= 1) { 517e5c89e4eSSatish Balay if (right == 1) { 518e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP2IntScalar(v[0],v[1],V[0],V[1],tmp,stmp); 519e5c89e4eSSatish Balay } 520e5c89e4eSSatish Balay PetscFunctionReturn(0); 521e5c89e4eSSatish Balay } 522e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[right/2],V[0],V[right/2],tmp,stmp); 523e5c89e4eSSatish Balay vl = v[0]; 524e5c89e4eSSatish Balay last = 0; 525e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 526e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP2IntScalar(v[last],v[i],V[last],V[i],tmp,stmp);} 527e5c89e4eSSatish Balay } 528e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[last],V[0],V[last],tmp,stmp); 529e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v,V,last-1);CHKERRQ(ierr); 530e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 531e5c89e4eSSatish Balay PetscFunctionReturn(0); 532e5c89e4eSSatish Balay } 533e5c89e4eSSatish Balay 534e5c89e4eSSatish Balay /*@ 535e5c89e4eSSatish Balay PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order; 536e5c89e4eSSatish Balay changes a second SCALAR array to match the sorted first INTEGER array. 537e5c89e4eSSatish Balay 538e5c89e4eSSatish Balay Not Collective 539e5c89e4eSSatish Balay 540e5c89e4eSSatish Balay Input Parameters: 541e5c89e4eSSatish Balay + n - number of values 542e5c89e4eSSatish Balay . i - array of integers 543e5c89e4eSSatish Balay - I - second array of scalars 544e5c89e4eSSatish Balay 545e5c89e4eSSatish Balay Level: intermediate 546e5c89e4eSSatish Balay 547e5c89e4eSSatish Balay Concepts: sorting^ints with array 548e5c89e4eSSatish Balay 549e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 550e5c89e4eSSatish Balay @*/ 5517087cfbeSBarry Smith PetscErrorCode PetscSortIntWithScalarArray(PetscInt n,PetscInt i[],PetscScalar Ii[]) 552e5c89e4eSSatish Balay { 553e5c89e4eSSatish Balay PetscErrorCode ierr; 554e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 555e5c89e4eSSatish Balay PetscScalar stmp; 556e5c89e4eSSatish Balay 557e5c89e4eSSatish Balay PetscFunctionBegin; 558e5c89e4eSSatish Balay if (n<8) { 559e5c89e4eSSatish Balay for (k=0; k<n; k++) { 560e5c89e4eSSatish Balay ik = i[k]; 561e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 562e5c89e4eSSatish Balay if (ik > i[j]) { 563b7940d39SSatish Balay SWAP2IntScalar(i[k],i[j],Ii[k],Ii[j],tmp,stmp); 564e5c89e4eSSatish Balay ik = i[k]; 565e5c89e4eSSatish Balay } 566e5c89e4eSSatish Balay } 567e5c89e4eSSatish Balay } 568e5c89e4eSSatish Balay } else { 569b7940d39SSatish Balay ierr = PetscSortIntWithScalarArray_Private(i,Ii,n-1);CHKERRQ(ierr); 570e5c89e4eSSatish Balay } 571e5c89e4eSSatish Balay PetscFunctionReturn(0); 572e5c89e4eSSatish Balay } 573e5c89e4eSSatish Balay 574ece8f864SToby Isaac #define SWAP2IntData(a,b,c,d,t,td,siz) \ 575ece8f864SToby Isaac do { \ 576ece8f864SToby Isaac PetscErrorCode _ierr; \ 577ece8f864SToby Isaac t=a;a=b;b=t; \ 578ece8f864SToby Isaac _ierr = PetscMemcpy(td,c,siz);CHKERRQ(_ierr); \ 579ece8f864SToby Isaac _ierr = PetscMemcpy(c,d,siz);CHKERRQ(_ierr); \ 580ece8f864SToby Isaac _ierr = PetscMemcpy(d,td,siz);CHKERRQ(_ierr); \ 581ece8f864SToby Isaac } while(0) 58217df18f2SToby Isaac 58317df18f2SToby Isaac /* 58417df18f2SToby Isaac Modified from PetscSortIntWithArray_Private(). 58517df18f2SToby Isaac */ 58617df18f2SToby Isaac static PetscErrorCode PetscSortIntWithDataArray_Private(PetscInt *v,char *V,PetscInt right,size_t size,void *work) 58717df18f2SToby Isaac { 58817df18f2SToby Isaac PetscErrorCode ierr; 58917df18f2SToby Isaac PetscInt i,vl,last,tmp; 59017df18f2SToby Isaac 59117df18f2SToby Isaac PetscFunctionBegin; 59217df18f2SToby Isaac if (right <= 1) { 59317df18f2SToby Isaac if (right == 1) { 59417df18f2SToby Isaac if (v[0] > v[1]) SWAP2IntData(v[0],v[1],V,V+size,tmp,work,size); 59517df18f2SToby Isaac } 59617df18f2SToby Isaac PetscFunctionReturn(0); 59717df18f2SToby Isaac } 59817df18f2SToby Isaac SWAP2IntData(v[0],v[right/2],V,V+size*(right/2),tmp,work,size); 59917df18f2SToby Isaac vl = v[0]; 60017df18f2SToby Isaac last = 0; 60117df18f2SToby Isaac for (i=1; i<=right; i++) { 60217df18f2SToby Isaac if (v[i] < vl) {last++; SWAP2IntData(v[last],v[i],V+size*last,V+size*i,tmp,work,size);} 60317df18f2SToby Isaac } 60417df18f2SToby Isaac SWAP2IntData(v[0],v[last],V,V + size*last,tmp,work,size); 60517df18f2SToby Isaac ierr = PetscSortIntWithDataArray_Private(v,V,last-1,size,work);CHKERRQ(ierr); 60617df18f2SToby Isaac ierr = PetscSortIntWithDataArray_Private(v+last+1,V+size*(last+1),right-(last+1),size,work);CHKERRQ(ierr); 60717df18f2SToby Isaac PetscFunctionReturn(0); 60817df18f2SToby Isaac } 60917df18f2SToby Isaac 6106c2863d0SBarry Smith /*@C 61117df18f2SToby Isaac PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order; 61217df18f2SToby Isaac changes a second array to match the sorted first INTEGER array. Unlike other sort routines, the user must 61317df18f2SToby Isaac provide workspace (the size of an element in the data array) to use when sorting. 61417df18f2SToby Isaac 61517df18f2SToby Isaac Not Collective 61617df18f2SToby Isaac 61717df18f2SToby Isaac Input Parameters: 61817df18f2SToby Isaac + n - number of values 61917df18f2SToby Isaac . i - array of integers 62017df18f2SToby Isaac . Ii - second array of data 62117df18f2SToby Isaac . size - sizeof elements in the data array in bytes 62217df18f2SToby Isaac - work - workspace of "size" bytes used when sorting 62317df18f2SToby Isaac 62417df18f2SToby Isaac Level: intermediate 62517df18f2SToby Isaac 62617df18f2SToby Isaac Concepts: sorting^ints with array 62717df18f2SToby Isaac 62817df18f2SToby Isaac .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 62917df18f2SToby Isaac @*/ 63017df18f2SToby Isaac PetscErrorCode PetscSortIntWithDataArray(PetscInt n,PetscInt i[],void *Ii,size_t size,void *work) 63117df18f2SToby Isaac { 63217df18f2SToby Isaac char *V = (char *) Ii; 63317df18f2SToby Isaac PetscErrorCode ierr; 63417df18f2SToby Isaac PetscInt j,k,tmp,ik; 63517df18f2SToby Isaac 63617df18f2SToby Isaac PetscFunctionBegin; 63717df18f2SToby Isaac if (n<8) { 63817df18f2SToby Isaac for (k=0; k<n; k++) { 63917df18f2SToby Isaac ik = i[k]; 64017df18f2SToby Isaac for (j=k+1; j<n; j++) { 64117df18f2SToby Isaac if (ik > i[j]) { 64217df18f2SToby Isaac SWAP2IntData(i[k],i[j],V+size*k,V+size*j,tmp,work,size); 64317df18f2SToby Isaac ik = i[k]; 64417df18f2SToby Isaac } 64517df18f2SToby Isaac } 64617df18f2SToby Isaac } 64717df18f2SToby Isaac } else { 64817df18f2SToby Isaac ierr = PetscSortIntWithDataArray_Private(i,V,n-1,size,work);CHKERRQ(ierr); 64917df18f2SToby Isaac } 65017df18f2SToby Isaac PetscFunctionReturn(0); 65117df18f2SToby Isaac } 65217df18f2SToby Isaac 653b4301105SBarry Smith 65421e72a00SBarry Smith /*@ 65521e72a00SBarry Smith PetscMergeIntArray - Merges two SORTED integer arrays, removes duplicate elements. 65621e72a00SBarry Smith 65721e72a00SBarry Smith Not Collective 65821e72a00SBarry Smith 65921e72a00SBarry Smith Input Parameters: 66021e72a00SBarry Smith + an - number of values in the first array 66121e72a00SBarry Smith . aI - first sorted array of integers 66221e72a00SBarry Smith . bn - number of values in the second array 66321e72a00SBarry Smith - bI - second array of integers 66421e72a00SBarry Smith 66521e72a00SBarry Smith Output Parameters: 66621e72a00SBarry Smith + n - number of values in the merged array 6676c2863d0SBarry Smith - L - merged sorted array, this is allocated if an array is not provided 66821e72a00SBarry Smith 66921e72a00SBarry Smith Level: intermediate 67021e72a00SBarry Smith 67121e72a00SBarry Smith Concepts: merging^arrays 67221e72a00SBarry Smith 67321e72a00SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 67421e72a00SBarry Smith @*/ 6756c2863d0SBarry Smith PetscErrorCode PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L) 67621e72a00SBarry Smith { 67721e72a00SBarry Smith PetscErrorCode ierr; 67821e72a00SBarry Smith PetscInt *L_ = *L, ak, bk, k; 67921e72a00SBarry Smith 68021e72a00SBarry Smith if (!L_) { 68121e72a00SBarry Smith ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr); 68221e72a00SBarry Smith L_ = *L; 68321e72a00SBarry Smith } 68421e72a00SBarry Smith k = ak = bk = 0; 68521e72a00SBarry Smith while (ak < an && bk < bn) { 68621e72a00SBarry Smith if (aI[ak] == bI[bk]) { 68721e72a00SBarry Smith L_[k] = aI[ak]; 68821e72a00SBarry Smith ++ak; 68921e72a00SBarry Smith ++bk; 69021e72a00SBarry Smith ++k; 69121e72a00SBarry Smith } else if (aI[ak] < bI[bk]) { 69221e72a00SBarry Smith L_[k] = aI[ak]; 69321e72a00SBarry Smith ++ak; 69421e72a00SBarry Smith ++k; 69521e72a00SBarry Smith } else { 69621e72a00SBarry Smith L_[k] = bI[bk]; 69721e72a00SBarry Smith ++bk; 69821e72a00SBarry Smith ++k; 69921e72a00SBarry Smith } 70021e72a00SBarry Smith } 70121e72a00SBarry Smith if (ak < an) { 70221e72a00SBarry Smith ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr); 70321e72a00SBarry Smith k += (an-ak); 70421e72a00SBarry Smith } 70521e72a00SBarry Smith if (bk < bn) { 70621e72a00SBarry Smith ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr); 70721e72a00SBarry Smith k += (bn-bk); 70821e72a00SBarry Smith } 70921e72a00SBarry Smith *n = k; 71021e72a00SBarry Smith PetscFunctionReturn(0); 71121e72a00SBarry Smith } 712b4301105SBarry Smith 713c1f0200aSDmitry Karpeev /*@ 71421e72a00SBarry Smith PetscMergeIntArrayPair - Merges two SORTED integer arrays that share NO common values along with an additional array of integers. 715c1f0200aSDmitry Karpeev The additional arrays are the same length as sorted arrays and are merged 716c1f0200aSDmitry Karpeev in the order determined by the merging of the sorted pair. 717c1f0200aSDmitry Karpeev 718c1f0200aSDmitry Karpeev Not Collective 719c1f0200aSDmitry Karpeev 720c1f0200aSDmitry Karpeev Input Parameters: 721c1f0200aSDmitry Karpeev + an - number of values in the first array 722c1f0200aSDmitry Karpeev . aI - first sorted array of integers 723c1f0200aSDmitry Karpeev . aJ - first additional array of integers 724c1f0200aSDmitry Karpeev . bn - number of values in the second array 725c1f0200aSDmitry Karpeev . bI - second array of integers 726c1f0200aSDmitry Karpeev - bJ - second additional of integers 727c1f0200aSDmitry Karpeev 728c1f0200aSDmitry Karpeev Output Parameters: 729c1f0200aSDmitry Karpeev + n - number of values in the merged array (== an + bn) 73014c49006SJed Brown . L - merged sorted array 731c1f0200aSDmitry Karpeev - J - merged additional array 732c1f0200aSDmitry Karpeev 73370d8d27cSBarry 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 734c1f0200aSDmitry Karpeev Level: intermediate 735c1f0200aSDmitry Karpeev 736c1f0200aSDmitry Karpeev Concepts: merging^arrays 737c1f0200aSDmitry Karpeev 738c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 739c1f0200aSDmitry Karpeev @*/ 7406c2863d0SBarry 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) 741c1f0200aSDmitry Karpeev { 742c1f0200aSDmitry Karpeev PetscErrorCode ierr; 74370d8d27cSBarry Smith PetscInt n_, *L_, *J_, ak, bk, k; 744c1f0200aSDmitry Karpeev 74514c49006SJed Brown PetscFunctionBegin; 74670d8d27cSBarry Smith PetscValidIntPointer(L,8); 74770d8d27cSBarry Smith PetscValidIntPointer(J,9); 748c1f0200aSDmitry Karpeev n_ = an + bn; 749c1f0200aSDmitry Karpeev *n = n_; 75070d8d27cSBarry Smith if (!*L) { 751785e854fSJed Brown ierr = PetscMalloc1(n_, L);CHKERRQ(ierr); 75270d8d27cSBarry Smith } 753d7aa01a8SHong Zhang L_ = *L; 75470d8d27cSBarry Smith if (!*J) { 75570d8d27cSBarry Smith ierr = PetscMalloc1(n_, J);CHKERRQ(ierr); 756c1f0200aSDmitry Karpeev } 757c1f0200aSDmitry Karpeev J_ = *J; 758c1f0200aSDmitry Karpeev k = ak = bk = 0; 759c1f0200aSDmitry Karpeev while (ak < an && bk < bn) { 760c1f0200aSDmitry Karpeev if (aI[ak] <= bI[bk]) { 761d7aa01a8SHong Zhang L_[k] = aI[ak]; 762c1f0200aSDmitry Karpeev J_[k] = aJ[ak]; 763c1f0200aSDmitry Karpeev ++ak; 764c1f0200aSDmitry Karpeev ++k; 765a297a907SKarl Rupp } else { 766d7aa01a8SHong Zhang L_[k] = bI[bk]; 767c1f0200aSDmitry Karpeev J_[k] = bJ[bk]; 768c1f0200aSDmitry Karpeev ++bk; 769c1f0200aSDmitry Karpeev ++k; 770c1f0200aSDmitry Karpeev } 771c1f0200aSDmitry Karpeev } 772c1f0200aSDmitry Karpeev if (ak < an) { 773d7aa01a8SHong Zhang ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr); 774c1f0200aSDmitry Karpeev ierr = PetscMemcpy(J_+k,aJ+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr); 775c1f0200aSDmitry Karpeev k += (an-ak); 776c1f0200aSDmitry Karpeev } 777c1f0200aSDmitry Karpeev if (bk < bn) { 778d7aa01a8SHong Zhang ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr); 779c1f0200aSDmitry Karpeev ierr = PetscMemcpy(J_+k,bJ+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr); 780c1f0200aSDmitry Karpeev } 781c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 782c1f0200aSDmitry Karpeev } 783c1f0200aSDmitry Karpeev 784e498c390SJed Brown /*@ 785e498c390SJed Brown PetscMergeMPIIntArray - Merges two SORTED integer arrays. 786e498c390SJed Brown 787e498c390SJed Brown Not Collective 788e498c390SJed Brown 789e498c390SJed Brown Input Parameters: 790e498c390SJed Brown + an - number of values in the first array 791e498c390SJed Brown . aI - first sorted array of integers 792e498c390SJed Brown . bn - number of values in the second array 793e498c390SJed Brown - bI - second array of integers 794e498c390SJed Brown 795e498c390SJed Brown Output Parameters: 796e498c390SJed Brown + n - number of values in the merged array (<= an + bn) 797e498c390SJed Brown - L - merged sorted array, allocated if address of NULL pointer is passed 798e498c390SJed Brown 799e498c390SJed Brown Level: intermediate 800e498c390SJed Brown 801e498c390SJed Brown Concepts: merging^arrays 802e498c390SJed Brown 803e498c390SJed Brown .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 804e498c390SJed Brown @*/ 805e498c390SJed Brown PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L) 806e498c390SJed Brown { 807e498c390SJed Brown PetscErrorCode ierr; 808e498c390SJed Brown PetscInt ai,bi,k; 809e498c390SJed Brown 810e498c390SJed Brown PetscFunctionBegin; 811e498c390SJed Brown if (!*L) {ierr = PetscMalloc1((an+bn),L);CHKERRQ(ierr);} 812e498c390SJed Brown for (ai=0,bi=0,k=0; ai<an || bi<bn; ) { 813e498c390SJed Brown PetscInt t = -1; 814e498c390SJed Brown for ( ; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai]; 815e498c390SJed Brown for ( ; bi<bn && bI[bi] == t; bi++); 816e498c390SJed Brown for ( ; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi]; 817e498c390SJed Brown for ( ; ai<an && aI[ai] == t; ai++); 818e498c390SJed Brown } 819e498c390SJed Brown *n = k; 820e498c390SJed Brown PetscFunctionReturn(0); 821e498c390SJed Brown } 822e5c89e4eSSatish Balay 8236c2863d0SBarry Smith /*@C 82442eaa7f3SBarry Smith PetscProcessTree - Prepares tree data to be displayed graphically 82542eaa7f3SBarry Smith 82642eaa7f3SBarry Smith Not Collective 82742eaa7f3SBarry Smith 82842eaa7f3SBarry Smith Input Parameters: 82942eaa7f3SBarry Smith + n - number of values 83042eaa7f3SBarry Smith . mask - indicates those entries in the tree, location 0 is always masked 83142eaa7f3SBarry Smith - parentid - indicates the parent of each entry 83242eaa7f3SBarry Smith 83342eaa7f3SBarry Smith Output Parameters: 83442eaa7f3SBarry Smith + Nlevels - the number of levels 83542eaa7f3SBarry Smith . Level - for each node tells its level 83642eaa7f3SBarry Smith . Levelcnts - the number of nodes on each level 83742eaa7f3SBarry Smith . Idbylevel - a list of ids on each of the levels, first level followed by second etc 83842eaa7f3SBarry Smith - Column - for each id tells its column index 83942eaa7f3SBarry Smith 8406c2863d0SBarry Smith Level: developer 84142eaa7f3SBarry Smith 8426c2863d0SBarry Smith Notes: This code is not currently used 84342eaa7f3SBarry Smith 84442eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation() 84542eaa7f3SBarry Smith @*/ 8467087cfbeSBarry Smith PetscErrorCode PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column) 84742eaa7f3SBarry Smith { 84842eaa7f3SBarry Smith PetscInt i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column; 84942eaa7f3SBarry Smith PetscErrorCode ierr; 850ace3abfcSBarry Smith PetscBool done = PETSC_FALSE; 85142eaa7f3SBarry Smith 85242eaa7f3SBarry Smith PetscFunctionBegin; 85342eaa7f3SBarry Smith if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set"); 85442eaa7f3SBarry Smith for (i=0; i<n; i++) { 85542eaa7f3SBarry Smith if (mask[i]) continue; 85642eaa7f3SBarry Smith if (parentid[i] == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent"); 85742eaa7f3SBarry Smith if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked"); 85842eaa7f3SBarry Smith } 85942eaa7f3SBarry Smith 86042eaa7f3SBarry Smith for (i=0; i<n; i++) { 86142eaa7f3SBarry Smith if (!mask[i]) nmask++; 86242eaa7f3SBarry Smith } 86342eaa7f3SBarry Smith 86442eaa7f3SBarry Smith /* determine the level in the tree of each node */ 8651795a4d1SJed Brown ierr = PetscCalloc1(n,&level);CHKERRQ(ierr); 866a297a907SKarl Rupp 86742eaa7f3SBarry Smith level[0] = 1; 86842eaa7f3SBarry Smith while (!done) { 86942eaa7f3SBarry Smith done = PETSC_TRUE; 87042eaa7f3SBarry Smith for (i=0; i<n; i++) { 87142eaa7f3SBarry Smith if (mask[i]) continue; 87242eaa7f3SBarry Smith if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1; 87342eaa7f3SBarry Smith else if (!level[i]) done = PETSC_FALSE; 87442eaa7f3SBarry Smith } 87542eaa7f3SBarry Smith } 87642eaa7f3SBarry Smith for (i=0; i<n; i++) { 87742eaa7f3SBarry Smith level[i]--; 87842eaa7f3SBarry Smith nlevels = PetscMax(nlevels,level[i]); 87942eaa7f3SBarry Smith } 88042eaa7f3SBarry Smith 88142eaa7f3SBarry Smith /* count the number of nodes on each level and its max */ 8821795a4d1SJed Brown ierr = PetscCalloc1(nlevels,&levelcnt);CHKERRQ(ierr); 88342eaa7f3SBarry Smith for (i=0; i<n; i++) { 88442eaa7f3SBarry Smith if (mask[i]) continue; 88542eaa7f3SBarry Smith levelcnt[level[i]-1]++; 88642eaa7f3SBarry Smith } 887a297a907SKarl Rupp for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]); 88842eaa7f3SBarry Smith 88942eaa7f3SBarry Smith /* for each level sort the ids by the parent id */ 890dcca6d9dSJed Brown ierr = PetscMalloc2(levelmax,&workid,levelmax,&workparentid);CHKERRQ(ierr); 891785e854fSJed Brown ierr = PetscMalloc1(nmask,&idbylevel);CHKERRQ(ierr); 89242eaa7f3SBarry Smith for (j=1; j<=nlevels;j++) { 89342eaa7f3SBarry Smith cnt = 0; 89442eaa7f3SBarry Smith for (i=0; i<n; i++) { 89542eaa7f3SBarry Smith if (mask[i]) continue; 89642eaa7f3SBarry Smith if (level[i] != j) continue; 89742eaa7f3SBarry Smith workid[cnt] = i; 89842eaa7f3SBarry Smith workparentid[cnt++] = parentid[i]; 89942eaa7f3SBarry Smith } 90042eaa7f3SBarry Smith /* PetscIntView(cnt,workparentid,0); 90142eaa7f3SBarry Smith PetscIntView(cnt,workid,0); 90242eaa7f3SBarry Smith ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr); 90342eaa7f3SBarry Smith PetscIntView(cnt,workparentid,0); 90442eaa7f3SBarry Smith PetscIntView(cnt,workid,0);*/ 90542eaa7f3SBarry Smith ierr = PetscMemcpy(idbylevel+tcnt,workid,cnt*sizeof(PetscInt));CHKERRQ(ierr); 90642eaa7f3SBarry Smith tcnt += cnt; 90742eaa7f3SBarry Smith } 90842eaa7f3SBarry Smith if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes"); 90942eaa7f3SBarry Smith ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr); 91042eaa7f3SBarry Smith 91142eaa7f3SBarry Smith /* for each node list its column */ 912785e854fSJed Brown ierr = PetscMalloc1(n,&column);CHKERRQ(ierr); 91342eaa7f3SBarry Smith cnt = 0; 91442eaa7f3SBarry Smith for (j=0; j<nlevels; j++) { 91542eaa7f3SBarry Smith for (i=0; i<levelcnt[j]; i++) { 91642eaa7f3SBarry Smith column[idbylevel[cnt++]] = i; 91742eaa7f3SBarry Smith } 91842eaa7f3SBarry Smith } 91942eaa7f3SBarry Smith 92042eaa7f3SBarry Smith *Nlevels = nlevels; 92142eaa7f3SBarry Smith *Level = level; 92242eaa7f3SBarry Smith *Levelcnt = levelcnt; 92342eaa7f3SBarry Smith *Idbylevel = idbylevel; 92442eaa7f3SBarry Smith *Column = column; 92542eaa7f3SBarry Smith PetscFunctionReturn(0); 92642eaa7f3SBarry Smith } 927