17d0a6c19SBarry Smith 2e5c89e4eSSatish Balay /* 3e5c89e4eSSatish Balay This file contains routines for sorting integers. Values are sorted in place. 4e5c89e4eSSatish Balay */ 5c6db04a5SJed Brown #include <petscsys.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 #undef __FUNCT__ 23e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortInt_Private" 24e5c89e4eSSatish Balay /* 25e5c89e4eSSatish Balay A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 26e5c89e4eSSatish Balay Assumes 0 origin for v, number of elements = right+1 (right is index of 27e5c89e4eSSatish Balay right-most member). 28e5c89e4eSSatish Balay */ 29ef8e3583SJed Brown static void PetscSortInt_Private(PetscInt *v,PetscInt right) 30e5c89e4eSSatish Balay { 31ef8e3583SJed Brown PetscInt i,j,pivot,tmp; 32e5c89e4eSSatish Balay 33e5c89e4eSSatish Balay if (right <= 1) { 34e5c89e4eSSatish Balay if (right == 1) { 35e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP(v[0],v[1],tmp); 36e5c89e4eSSatish Balay } 37ef8e3583SJed Brown return; 38e5c89e4eSSatish Balay } 39ef8e3583SJed Brown i = MEDIAN(v,right); /* Choose a pivot */ 40ef8e3583SJed Brown SWAP(v[0],v[i],tmp); /* Move it out of the way */ 41ef8e3583SJed Brown pivot = v[0]; 42ef8e3583SJed Brown for (i=0,j=right+1;;) { 43ef8e3583SJed Brown while (++i < j && v[i] <= pivot); /* Scan from the left */ 44ef8e3583SJed Brown while (v[--j] > pivot); /* Scan from the right */ 45ef8e3583SJed Brown if (i >= j) break; 46ef8e3583SJed Brown SWAP(v[i],v[j],tmp); 47e5c89e4eSSatish Balay } 48ef8e3583SJed Brown SWAP(v[0],v[j],tmp); /* Put pivot back in place. */ 49ef8e3583SJed Brown PetscSortInt_Private(v,j-1); 50ef8e3583SJed Brown PetscSortInt_Private(v+j+1,right-(j+1)); 51e5c89e4eSSatish Balay } 52e5c89e4eSSatish Balay 53e5c89e4eSSatish Balay #undef __FUNCT__ 54e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortInt" 55e5c89e4eSSatish Balay /*@ 56e5c89e4eSSatish Balay PetscSortInt - Sorts an array of integers in place in increasing order. 57e5c89e4eSSatish Balay 58e5c89e4eSSatish Balay Not Collective 59e5c89e4eSSatish Balay 60e5c89e4eSSatish Balay Input Parameters: 61e5c89e4eSSatish Balay + n - number of values 62e5c89e4eSSatish Balay - i - array of integers 63e5c89e4eSSatish Balay 64e5c89e4eSSatish Balay Level: intermediate 65e5c89e4eSSatish Balay 66e5c89e4eSSatish Balay Concepts: sorting^ints 67e5c89e4eSSatish Balay 68e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntWithPermutation() 69e5c89e4eSSatish Balay @*/ 707087cfbeSBarry Smith PetscErrorCode PetscSortInt(PetscInt n,PetscInt i[]) 71e5c89e4eSSatish Balay { 72e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 73e5c89e4eSSatish Balay 74e5c89e4eSSatish Balay PetscFunctionBegin; 75e5c89e4eSSatish Balay if (n<8) { 76e5c89e4eSSatish Balay for (k=0; k<n; k++) { 77e5c89e4eSSatish Balay ik = i[k]; 78e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 79e5c89e4eSSatish Balay if (ik > i[j]) { 80e5c89e4eSSatish Balay SWAP(i[k],i[j],tmp); 81e5c89e4eSSatish Balay ik = i[k]; 82e5c89e4eSSatish Balay } 83e5c89e4eSSatish Balay } 84e5c89e4eSSatish Balay } 85e5c89e4eSSatish Balay } else { 86ef8e3583SJed Brown PetscSortInt_Private(i,n-1); 87e5c89e4eSSatish Balay } 88e5c89e4eSSatish Balay PetscFunctionReturn(0); 89e5c89e4eSSatish Balay } 90e5c89e4eSSatish Balay 911db36a52SBarry Smith #undef __FUNCT__ 921db36a52SBarry Smith #define __FUNCT__ "PetscSortRemoveDupsInt" 931db36a52SBarry Smith /*@ 941db36a52SBarry Smith PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries 951db36a52SBarry Smith 961db36a52SBarry Smith Not Collective 971db36a52SBarry Smith 981db36a52SBarry Smith Input Parameters: 991db36a52SBarry Smith + n - number of values 1001db36a52SBarry Smith - ii - array of integers 1011db36a52SBarry Smith 1021db36a52SBarry Smith Output Parameter: 1031db36a52SBarry Smith . n - number of non-redundant values 1041db36a52SBarry Smith 1051db36a52SBarry Smith Level: intermediate 1061db36a52SBarry Smith 1071db36a52SBarry Smith Concepts: sorting^ints 1081db36a52SBarry Smith 1091db36a52SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt() 1101db36a52SBarry Smith @*/ 1117087cfbeSBarry Smith PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n,PetscInt ii[]) 1121db36a52SBarry Smith { 1131db36a52SBarry Smith PetscErrorCode ierr; 1141db36a52SBarry Smith PetscInt i,s = 0,N = *n, b = 0; 1151db36a52SBarry Smith 1161db36a52SBarry Smith PetscFunctionBegin; 1171db36a52SBarry Smith ierr = PetscSortInt(N,ii);CHKERRQ(ierr); 1181db36a52SBarry Smith for (i=0; i<N-1; i++) { 119a5891931SBarry Smith if (ii[b+s+1] != ii[b]) {ii[b+1] = ii[b+s+1]; b++;} 120a5891931SBarry Smith else s++; 1211db36a52SBarry Smith } 1221db36a52SBarry Smith *n = N - s; 1231db36a52SBarry Smith PetscFunctionReturn(0); 1241db36a52SBarry Smith } 1251db36a52SBarry Smith 12660e03357SMatthew G Knepley #undef __FUNCT__ 12760e03357SMatthew G Knepley #define __FUNCT__ "PetscFindInt" 12860e03357SMatthew G Knepley /*@ 12960e03357SMatthew G Knepley PetscFindInt - Finds and integers in a sorted array of integers 13060e03357SMatthew G Knepley 13160e03357SMatthew G Knepley Not Collective 13260e03357SMatthew G Knepley 13360e03357SMatthew G Knepley Input Parameters: 13460e03357SMatthew G Knepley + key - the integer to locate 13560e03357SMatthew G Knepley . n - number of value in the array 13660e03357SMatthew G Knepley - ii - array of integers 13760e03357SMatthew G Knepley 13860e03357SMatthew G Knepley Output Parameter: 13960e03357SMatthew G Knepley . loc - the location, or -1 if the value is not found 14060e03357SMatthew G Knepley 14160e03357SMatthew G Knepley Level: intermediate 14260e03357SMatthew G Knepley 14360e03357SMatthew G Knepley Concepts: sorting^ints 14460e03357SMatthew G Knepley 14560e03357SMatthew G Knepley .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt() 14660e03357SMatthew G Knepley @*/ 147*026ec6cbSMatthew G Knepley PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt ii[], PetscInt *loc) 14860e03357SMatthew G Knepley { 14960e03357SMatthew G Knepley /* inclusive indices 15060e03357SMatthew G Knepley 0 <= imin when using truncate toward zero divide 15160e03357SMatthew G Knepley imid = (imin+imax)/2; 15260e03357SMatthew G Knepley imin unrestricted when using truncate toward minus infinity divide 15360e03357SMatthew G Knepley imid = (imin+imax)>>1; or 15460e03357SMatthew G Knepley imid = (int)floor((imin+imax)/2.0); */ 15560e03357SMatthew G Knepley PetscInt imin = 0, imax = n-1; 15660e03357SMatthew G Knepley 15760e03357SMatthew G Knepley PetscFunctionBegin; 15860e03357SMatthew G Knepley PetscValidPointer(ii, 3); 15960e03357SMatthew G Knepley PetscValidPointer(loc, 4); 16060e03357SMatthew G Knepley /* continually narrow search until just one element remains */ 16160e03357SMatthew G Knepley while(imin < imax) { 16260e03357SMatthew G Knepley PetscInt imid = (imin+imax)/2; 16360e03357SMatthew G Knepley 16460e03357SMatthew G Knepley /* code must guarantee the interval is reduced at each iteration 16560e03357SMatthew G Knepley assert(imid < imax); */ 16660e03357SMatthew G Knepley /* reduce the search */ 16760e03357SMatthew G Knepley if (ii[imid] < key) { 16860e03357SMatthew G Knepley imin = imid + 1; 16960e03357SMatthew G Knepley } else { 17060e03357SMatthew G Knepley imax = imid; 17160e03357SMatthew G Knepley } 17260e03357SMatthew G Knepley } 17360e03357SMatthew G Knepley /* At exit of while: 17460e03357SMatthew G Knepley if ii[] is empty, then imax < imin 17560e03357SMatthew G Knepley otherwise imax == imin */ 17660e03357SMatthew G Knepley /* deferred test for equality */ 17760e03357SMatthew G Knepley if ((imax == imin) && (ii[imin] == key)) { 17860e03357SMatthew G Knepley *loc = imin; 17960e03357SMatthew G Knepley } else { 18060e03357SMatthew G Knepley *loc = -1; 18160e03357SMatthew G Knepley } 18260e03357SMatthew G Knepley PetscFunctionReturn(0); 18360e03357SMatthew G Knepley } 18460e03357SMatthew G Knepley 1851db36a52SBarry Smith 186e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 187e5c89e4eSSatish Balay #define SWAP2(a,b,c,d,t) {t=a;a=b;b=t;t=c;c=d;d=t;} 188e5c89e4eSSatish Balay 189e5c89e4eSSatish Balay #undef __FUNCT__ 190e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithArray_Private" 191e5c89e4eSSatish Balay /* 192e5c89e4eSSatish Balay A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 193e5c89e4eSSatish Balay Assumes 0 origin for v, number of elements = right+1 (right is index of 194e5c89e4eSSatish Balay right-most member). 195e5c89e4eSSatish Balay */ 196e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithArray_Private(PetscInt *v,PetscInt *V,PetscInt right) 197e5c89e4eSSatish Balay { 198e5c89e4eSSatish Balay PetscErrorCode ierr; 199e5c89e4eSSatish Balay PetscInt i,vl,last,tmp; 200e5c89e4eSSatish Balay 201e5c89e4eSSatish Balay PetscFunctionBegin; 202e5c89e4eSSatish Balay if (right <= 1) { 203e5c89e4eSSatish Balay if (right == 1) { 204e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp); 205e5c89e4eSSatish Balay } 206e5c89e4eSSatish Balay PetscFunctionReturn(0); 207e5c89e4eSSatish Balay } 208e5c89e4eSSatish Balay SWAP2(v[0],v[right/2],V[0],V[right/2],tmp); 209e5c89e4eSSatish Balay vl = v[0]; 210e5c89e4eSSatish Balay last = 0; 211e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 212e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);} 213e5c89e4eSSatish Balay } 214e5c89e4eSSatish Balay SWAP2(v[0],v[last],V[0],V[last],tmp); 215e5c89e4eSSatish Balay ierr = PetscSortIntWithArray_Private(v,V,last-1);CHKERRQ(ierr); 216e5c89e4eSSatish Balay ierr = PetscSortIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 217e5c89e4eSSatish Balay PetscFunctionReturn(0); 218e5c89e4eSSatish Balay } 219e5c89e4eSSatish Balay 220e5c89e4eSSatish Balay #undef __FUNCT__ 221e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithArray" 222e5c89e4eSSatish Balay /*@ 223e5c89e4eSSatish Balay PetscSortIntWithArray - Sorts an array of integers in place in increasing order; 224e5c89e4eSSatish Balay changes a second array to match the sorted first array. 225e5c89e4eSSatish Balay 226e5c89e4eSSatish Balay Not Collective 227e5c89e4eSSatish Balay 228e5c89e4eSSatish Balay Input Parameters: 229e5c89e4eSSatish Balay + n - number of values 230e5c89e4eSSatish Balay . i - array of integers 231e5c89e4eSSatish Balay - I - second array of integers 232e5c89e4eSSatish Balay 233e5c89e4eSSatish Balay Level: intermediate 234e5c89e4eSSatish Balay 235e5c89e4eSSatish Balay Concepts: sorting^ints with array 236e5c89e4eSSatish Balay 237e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 238e5c89e4eSSatish Balay @*/ 2397087cfbeSBarry Smith PetscErrorCode PetscSortIntWithArray(PetscInt n,PetscInt i[],PetscInt Ii[]) 240e5c89e4eSSatish Balay { 241e5c89e4eSSatish Balay PetscErrorCode ierr; 242e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 243e5c89e4eSSatish Balay 244e5c89e4eSSatish Balay PetscFunctionBegin; 245e5c89e4eSSatish Balay if (n<8) { 246e5c89e4eSSatish Balay for (k=0; k<n; k++) { 247e5c89e4eSSatish Balay ik = i[k]; 248e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 249e5c89e4eSSatish Balay if (ik > i[j]) { 250b7940d39SSatish Balay SWAP2(i[k],i[j],Ii[k],Ii[j],tmp); 251e5c89e4eSSatish Balay ik = i[k]; 252e5c89e4eSSatish Balay } 253e5c89e4eSSatish Balay } 254e5c89e4eSSatish Balay } 255e5c89e4eSSatish Balay } else { 256b7940d39SSatish Balay ierr = PetscSortIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr); 257e5c89e4eSSatish Balay } 258e5c89e4eSSatish Balay PetscFunctionReturn(0); 259e5c89e4eSSatish Balay } 260e5c89e4eSSatish Balay 261c1f0200aSDmitry Karpeev 262c1f0200aSDmitry 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;} 263c1f0200aSDmitry Karpeev 264c1f0200aSDmitry Karpeev #undef __FUNCT__ 265c1f0200aSDmitry Karpeev #define __FUNCT__ "PetscSortIntWithArrayPair_Private" 266c1f0200aSDmitry Karpeev /* 267c1f0200aSDmitry Karpeev A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 268c1f0200aSDmitry Karpeev Assumes 0 origin for v, number of elements = right+1 (right is index of 269c1f0200aSDmitry Karpeev right-most member). 270c1f0200aSDmitry Karpeev */ 271d7aa01a8SHong Zhang static PetscErrorCode PetscSortIntWithArrayPair_Private(PetscInt *L,PetscInt *J, PetscInt *K, PetscInt right) 272c1f0200aSDmitry Karpeev { 273c1f0200aSDmitry Karpeev PetscErrorCode ierr; 274c1f0200aSDmitry Karpeev PetscInt i,vl,last,tmp; 275c1f0200aSDmitry Karpeev 276c1f0200aSDmitry Karpeev PetscFunctionBegin; 277c1f0200aSDmitry Karpeev if (right <= 1) { 278c1f0200aSDmitry Karpeev if (right == 1) { 279d7aa01a8SHong Zhang if (L[0] > L[1]) SWAP3(L[0],L[1],J[0],J[1],K[0],K[1],tmp); 280c1f0200aSDmitry Karpeev } 281c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 282c1f0200aSDmitry Karpeev } 283d7aa01a8SHong Zhang SWAP3(L[0],L[right/2],J[0],J[right/2],K[0],K[right/2],tmp); 284d7aa01a8SHong Zhang vl = L[0]; 285c1f0200aSDmitry Karpeev last = 0; 286c1f0200aSDmitry Karpeev for (i=1; i<=right; i++) { 287d7aa01a8SHong Zhang if (L[i] < vl) {last++; SWAP3(L[last],L[i],J[last],J[i],K[last],K[i],tmp);} 288c1f0200aSDmitry Karpeev } 289d7aa01a8SHong Zhang SWAP3(L[0],L[last],J[0],J[last],K[0],K[last],tmp); 290d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L,J,K,last-1);CHKERRQ(ierr); 291d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L+last+1,J+last+1,K+last+1,right-(last+1));CHKERRQ(ierr); 292c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 293c1f0200aSDmitry Karpeev } 294c1f0200aSDmitry Karpeev 295c1f0200aSDmitry Karpeev #undef __FUNCT__ 296c1f0200aSDmitry Karpeev #define __FUNCT__ "PetscSortIntWithArrayPair" 297c1f0200aSDmitry Karpeev /*@ 298c1f0200aSDmitry Karpeev PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order; 299c1f0200aSDmitry Karpeev changes a pair of integer arrays to match the sorted first array. 300c1f0200aSDmitry Karpeev 301c1f0200aSDmitry Karpeev Not Collective 302c1f0200aSDmitry Karpeev 303c1f0200aSDmitry Karpeev Input Parameters: 304c1f0200aSDmitry Karpeev + n - number of values 305c1f0200aSDmitry Karpeev . I - array of integers 306c1f0200aSDmitry Karpeev . J - second array of integers (first array of the pair) 307c1f0200aSDmitry Karpeev - K - third array of integers (second array of the pair) 308c1f0200aSDmitry Karpeev 309c1f0200aSDmitry Karpeev Level: intermediate 310c1f0200aSDmitry Karpeev 311c1f0200aSDmitry Karpeev Concepts: sorting^ints with array pair 312c1f0200aSDmitry Karpeev 313c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray() 314c1f0200aSDmitry Karpeev @*/ 315d7aa01a8SHong Zhang PetscErrorCode PetscSortIntWithArrayPair(PetscInt n,PetscInt *L,PetscInt *J, PetscInt *K) 316c1f0200aSDmitry Karpeev { 317c1f0200aSDmitry Karpeev PetscErrorCode ierr; 318c1f0200aSDmitry Karpeev PetscInt j,k,tmp,ik; 319c1f0200aSDmitry Karpeev 320c1f0200aSDmitry Karpeev PetscFunctionBegin; 321c1f0200aSDmitry Karpeev if (n<8) { 322c1f0200aSDmitry Karpeev for (k=0; k<n; k++) { 323d7aa01a8SHong Zhang ik = L[k]; 324c1f0200aSDmitry Karpeev for (j=k+1; j<n; j++) { 325d7aa01a8SHong Zhang if (ik > L[j]) { 326d7aa01a8SHong Zhang SWAP3(L[k],L[j],J[k],J[j],K[k],K[j],tmp); 327d7aa01a8SHong Zhang ik = L[k]; 328c1f0200aSDmitry Karpeev } 329c1f0200aSDmitry Karpeev } 330c1f0200aSDmitry Karpeev } 331c1f0200aSDmitry Karpeev } else { 332d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L,J,K,n-1);CHKERRQ(ierr); 333c1f0200aSDmitry Karpeev } 334c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 335c1f0200aSDmitry Karpeev } 336c1f0200aSDmitry Karpeev 33717d7d925SStefano Zampini #undef __FUNCT__ 33817d7d925SStefano Zampini #define __FUNCT__ "PetscSortMPIInt_Private" 33917d7d925SStefano Zampini /* 34017d7d925SStefano Zampini A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 34117d7d925SStefano Zampini Assumes 0 origin for v, number of elements = right+1 (right is index of 34217d7d925SStefano Zampini right-most member). 34317d7d925SStefano Zampini */ 34417d7d925SStefano Zampini static void PetscSortMPIInt_Private(PetscMPIInt *v,PetscInt right) 34517d7d925SStefano Zampini { 34617d7d925SStefano Zampini PetscInt i,j; 34717d7d925SStefano Zampini PetscMPIInt pivot,tmp; 34817d7d925SStefano Zampini 34917d7d925SStefano Zampini if (right <= 1) { 35017d7d925SStefano Zampini if (right == 1) { 35117d7d925SStefano Zampini if (v[0] > v[1]) SWAP(v[0],v[1],tmp); 35217d7d925SStefano Zampini } 35317d7d925SStefano Zampini return; 35417d7d925SStefano Zampini } 35517d7d925SStefano Zampini i = MEDIAN(v,right); /* Choose a pivot */ 35617d7d925SStefano Zampini SWAP(v[0],v[i],tmp); /* Move it out of the way */ 35717d7d925SStefano Zampini pivot = v[0]; 35817d7d925SStefano Zampini for (i=0,j=right+1;;) { 35917d7d925SStefano Zampini while (++i < j && v[i] <= pivot); /* Scan from the left */ 36017d7d925SStefano Zampini while (v[--j] > pivot); /* Scan from the right */ 36117d7d925SStefano Zampini if (i >= j) break; 36217d7d925SStefano Zampini SWAP(v[i],v[j],tmp); 36317d7d925SStefano Zampini } 36417d7d925SStefano Zampini SWAP(v[0],v[j],tmp); /* Put pivot back in place. */ 36517d7d925SStefano Zampini PetscSortMPIInt_Private(v,j-1); 36617d7d925SStefano Zampini PetscSortMPIInt_Private(v+j+1,right-(j+1)); 36717d7d925SStefano Zampini } 36817d7d925SStefano Zampini 36917d7d925SStefano Zampini #undef __FUNCT__ 37017d7d925SStefano Zampini #define __FUNCT__ "PetscSortMPIInt" 37117d7d925SStefano Zampini /*@ 37217d7d925SStefano Zampini PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order. 37317d7d925SStefano Zampini 37417d7d925SStefano Zampini Not Collective 37517d7d925SStefano Zampini 37617d7d925SStefano Zampini Input Parameters: 37717d7d925SStefano Zampini + n - number of values 37817d7d925SStefano Zampini - i - array of integers 37917d7d925SStefano Zampini 38017d7d925SStefano Zampini Level: intermediate 38117d7d925SStefano Zampini 38217d7d925SStefano Zampini Concepts: sorting^ints 38317d7d925SStefano Zampini 38417d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation() 38517d7d925SStefano Zampini @*/ 38617d7d925SStefano Zampini PetscErrorCode PetscSortMPIInt(PetscInt n,PetscMPIInt i[]) 38717d7d925SStefano Zampini { 38817d7d925SStefano Zampini PetscInt j,k; 38917d7d925SStefano Zampini PetscMPIInt tmp,ik; 39017d7d925SStefano Zampini 39117d7d925SStefano Zampini PetscFunctionBegin; 39217d7d925SStefano Zampini if (n<8) { 39317d7d925SStefano Zampini for (k=0; k<n; k++) { 39417d7d925SStefano Zampini ik = i[k]; 39517d7d925SStefano Zampini for (j=k+1; j<n; j++) { 39617d7d925SStefano Zampini if (ik > i[j]) { 39717d7d925SStefano Zampini SWAP(i[k],i[j],tmp); 39817d7d925SStefano Zampini ik = i[k]; 39917d7d925SStefano Zampini } 40017d7d925SStefano Zampini } 40117d7d925SStefano Zampini } 40217d7d925SStefano Zampini } else { 40317d7d925SStefano Zampini PetscSortMPIInt_Private(i,n-1); 40417d7d925SStefano Zampini } 40517d7d925SStefano Zampini PetscFunctionReturn(0); 40617d7d925SStefano Zampini } 40717d7d925SStefano Zampini 40817d7d925SStefano Zampini #undef __FUNCT__ 40917d7d925SStefano Zampini #define __FUNCT__ "PetscSortRemoveDupsMPIInt" 41017d7d925SStefano Zampini /*@ 41117d7d925SStefano Zampini PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries 41217d7d925SStefano Zampini 41317d7d925SStefano Zampini Not Collective 41417d7d925SStefano Zampini 41517d7d925SStefano Zampini Input Parameters: 41617d7d925SStefano Zampini + n - number of values 41717d7d925SStefano Zampini - ii - array of integers 41817d7d925SStefano Zampini 41917d7d925SStefano Zampini Output Parameter: 42017d7d925SStefano Zampini . n - number of non-redundant values 42117d7d925SStefano Zampini 42217d7d925SStefano Zampini Level: intermediate 42317d7d925SStefano Zampini 42417d7d925SStefano Zampini Concepts: sorting^ints 42517d7d925SStefano Zampini 42617d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt() 42717d7d925SStefano Zampini @*/ 42817d7d925SStefano Zampini PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt ii[]) 42917d7d925SStefano Zampini { 43017d7d925SStefano Zampini PetscErrorCode ierr; 43117d7d925SStefano Zampini PetscInt i,s = 0,N = *n, b = 0; 43217d7d925SStefano Zampini 43317d7d925SStefano Zampini PetscFunctionBegin; 43417d7d925SStefano Zampini ierr = PetscSortMPIInt(N,ii);CHKERRQ(ierr); 43517d7d925SStefano Zampini for (i=0; i<N-1; i++) { 43617d7d925SStefano Zampini if (ii[b+s+1] != ii[b]) {ii[b+1] = ii[b+s+1]; b++;} 43717d7d925SStefano Zampini else s++; 43817d7d925SStefano Zampini } 43917d7d925SStefano Zampini *n = N - s; 44017d7d925SStefano Zampini PetscFunctionReturn(0); 44117d7d925SStefano Zampini } 442c1f0200aSDmitry Karpeev 4434d615ea0SBarry Smith #undef __FUNCT__ 4444d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray_Private" 4454d615ea0SBarry Smith /* 4464d615ea0SBarry Smith A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 4474d615ea0SBarry Smith Assumes 0 origin for v, number of elements = right+1 (right is index of 4484d615ea0SBarry Smith right-most member). 4494d615ea0SBarry Smith */ 4504d615ea0SBarry Smith static PetscErrorCode PetscSortMPIIntWithArray_Private(PetscMPIInt *v,PetscMPIInt *V,PetscMPIInt right) 4514d615ea0SBarry Smith { 4524d615ea0SBarry Smith PetscErrorCode ierr; 4534d615ea0SBarry Smith PetscMPIInt i,vl,last,tmp; 4544d615ea0SBarry Smith 4554d615ea0SBarry Smith PetscFunctionBegin; 4564d615ea0SBarry Smith if (right <= 1) { 4574d615ea0SBarry Smith if (right == 1) { 4584d615ea0SBarry Smith if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp); 4594d615ea0SBarry Smith } 4604d615ea0SBarry Smith PetscFunctionReturn(0); 4614d615ea0SBarry Smith } 4624d615ea0SBarry Smith SWAP2(v[0],v[right/2],V[0],V[right/2],tmp); 4634d615ea0SBarry Smith vl = v[0]; 4644d615ea0SBarry Smith last = 0; 4654d615ea0SBarry Smith for (i=1; i<=right; i++) { 4664d615ea0SBarry Smith if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);} 4674d615ea0SBarry Smith } 4684d615ea0SBarry Smith SWAP2(v[0],v[last],V[0],V[last],tmp); 4694d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v,V,last-1);CHKERRQ(ierr); 4704d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 4714d615ea0SBarry Smith PetscFunctionReturn(0); 4724d615ea0SBarry Smith } 4734d615ea0SBarry Smith 4744d615ea0SBarry Smith #undef __FUNCT__ 4754d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray" 4764d615ea0SBarry Smith /*@ 4774d615ea0SBarry Smith PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order; 4784d615ea0SBarry Smith changes a second array to match the sorted first array. 4794d615ea0SBarry Smith 4804d615ea0SBarry Smith Not Collective 4814d615ea0SBarry Smith 4824d615ea0SBarry Smith Input Parameters: 4834d615ea0SBarry Smith + n - number of values 4844d615ea0SBarry Smith . i - array of integers 4854d615ea0SBarry Smith - I - second array of integers 4864d615ea0SBarry Smith 4874d615ea0SBarry Smith Level: intermediate 4884d615ea0SBarry Smith 4894d615ea0SBarry Smith Concepts: sorting^ints with array 4904d615ea0SBarry Smith 4914d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 4924d615ea0SBarry Smith @*/ 4937087cfbeSBarry Smith PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt i[],PetscMPIInt Ii[]) 4944d615ea0SBarry Smith { 4954d615ea0SBarry Smith PetscErrorCode ierr; 4964d615ea0SBarry Smith PetscMPIInt j,k,tmp,ik; 4974d615ea0SBarry Smith 4984d615ea0SBarry Smith PetscFunctionBegin; 4994d615ea0SBarry Smith if (n<8) { 5004d615ea0SBarry Smith for (k=0; k<n; k++) { 5014d615ea0SBarry Smith ik = i[k]; 5024d615ea0SBarry Smith for (j=k+1; j<n; j++) { 5034d615ea0SBarry Smith if (ik > i[j]) { 5044d615ea0SBarry Smith SWAP2(i[k],i[j],Ii[k],Ii[j],tmp); 5054d615ea0SBarry Smith ik = i[k]; 5064d615ea0SBarry Smith } 5074d615ea0SBarry Smith } 5084d615ea0SBarry Smith } 5094d615ea0SBarry Smith } else { 5104d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr); 5114d615ea0SBarry Smith } 5124d615ea0SBarry Smith PetscFunctionReturn(0); 5134d615ea0SBarry Smith } 5144d615ea0SBarry Smith 515e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 516e5c89e4eSSatish Balay #define SWAP2IntScalar(a,b,c,d,t,ts) {t=a;a=b;b=t;ts=c;c=d;d=ts;} 517e5c89e4eSSatish Balay 518e5c89e4eSSatish Balay #undef __FUNCT__ 519e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray_Private" 520e5c89e4eSSatish Balay /* 521e5c89e4eSSatish Balay Modified from PetscSortIntWithArray_Private(). 522e5c89e4eSSatish Balay */ 523e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithScalarArray_Private(PetscInt *v,PetscScalar *V,PetscInt right) 524e5c89e4eSSatish Balay { 525e5c89e4eSSatish Balay PetscErrorCode ierr; 526e5c89e4eSSatish Balay PetscInt i,vl,last,tmp; 527e5c89e4eSSatish Balay PetscScalar stmp; 528e5c89e4eSSatish Balay 529e5c89e4eSSatish Balay PetscFunctionBegin; 530e5c89e4eSSatish Balay if (right <= 1) { 531e5c89e4eSSatish Balay if (right == 1) { 532e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP2IntScalar(v[0],v[1],V[0],V[1],tmp,stmp); 533e5c89e4eSSatish Balay } 534e5c89e4eSSatish Balay PetscFunctionReturn(0); 535e5c89e4eSSatish Balay } 536e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[right/2],V[0],V[right/2],tmp,stmp); 537e5c89e4eSSatish Balay vl = v[0]; 538e5c89e4eSSatish Balay last = 0; 539e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 540e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP2IntScalar(v[last],v[i],V[last],V[i],tmp,stmp);} 541e5c89e4eSSatish Balay } 542e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[last],V[0],V[last],tmp,stmp); 543e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v,V,last-1);CHKERRQ(ierr); 544e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 545e5c89e4eSSatish Balay PetscFunctionReturn(0); 546e5c89e4eSSatish Balay } 547e5c89e4eSSatish Balay 548e5c89e4eSSatish Balay #undef __FUNCT__ 549e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray" 550e5c89e4eSSatish Balay /*@ 551e5c89e4eSSatish Balay PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order; 552e5c89e4eSSatish Balay changes a second SCALAR array to match the sorted first INTEGER array. 553e5c89e4eSSatish Balay 554e5c89e4eSSatish Balay Not Collective 555e5c89e4eSSatish Balay 556e5c89e4eSSatish Balay Input Parameters: 557e5c89e4eSSatish Balay + n - number of values 558e5c89e4eSSatish Balay . i - array of integers 559e5c89e4eSSatish Balay - I - second array of scalars 560e5c89e4eSSatish Balay 561e5c89e4eSSatish Balay Level: intermediate 562e5c89e4eSSatish Balay 563e5c89e4eSSatish Balay Concepts: sorting^ints with array 564e5c89e4eSSatish Balay 565e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 566e5c89e4eSSatish Balay @*/ 5677087cfbeSBarry Smith PetscErrorCode PetscSortIntWithScalarArray(PetscInt n,PetscInt i[],PetscScalar Ii[]) 568e5c89e4eSSatish Balay { 569e5c89e4eSSatish Balay PetscErrorCode ierr; 570e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 571e5c89e4eSSatish Balay PetscScalar stmp; 572e5c89e4eSSatish Balay 573e5c89e4eSSatish Balay PetscFunctionBegin; 574e5c89e4eSSatish Balay if (n<8) { 575e5c89e4eSSatish Balay for (k=0; k<n; k++) { 576e5c89e4eSSatish Balay ik = i[k]; 577e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 578e5c89e4eSSatish Balay if (ik > i[j]) { 579b7940d39SSatish Balay SWAP2IntScalar(i[k],i[j],Ii[k],Ii[j],tmp,stmp); 580e5c89e4eSSatish Balay ik = i[k]; 581e5c89e4eSSatish Balay } 582e5c89e4eSSatish Balay } 583e5c89e4eSSatish Balay } 584e5c89e4eSSatish Balay } else { 585b7940d39SSatish Balay ierr = PetscSortIntWithScalarArray_Private(i,Ii,n-1);CHKERRQ(ierr); 586e5c89e4eSSatish Balay } 587e5c89e4eSSatish Balay PetscFunctionReturn(0); 588e5c89e4eSSatish Balay } 589e5c89e4eSSatish Balay 590b4301105SBarry Smith 591b4301105SBarry Smith 592b4301105SBarry Smith #undef __FUNCT__ 593c1f0200aSDmitry Karpeev #define __FUNCT__ "PetscMergeIntArrayPair" 594c1f0200aSDmitry Karpeev /*@ 595c1f0200aSDmitry Karpeev PetscMergeIntArrayPair - Merges two SORTED integer arrays along with an additional array of integers. 596c1f0200aSDmitry Karpeev The additional arrays are the same length as sorted arrays and are merged 597c1f0200aSDmitry Karpeev in the order determined by the merging of the sorted pair. 598c1f0200aSDmitry Karpeev 599c1f0200aSDmitry Karpeev Not Collective 600c1f0200aSDmitry Karpeev 601c1f0200aSDmitry Karpeev Input Parameters: 602c1f0200aSDmitry Karpeev + an - number of values in the first array 603c1f0200aSDmitry Karpeev . aI - first sorted array of integers 604c1f0200aSDmitry Karpeev . aJ - first additional array of integers 605c1f0200aSDmitry Karpeev . bn - number of values in the second array 606c1f0200aSDmitry Karpeev . bI - second array of integers 607c1f0200aSDmitry Karpeev - bJ - second additional of integers 608c1f0200aSDmitry Karpeev 609c1f0200aSDmitry Karpeev Output Parameters: 610c1f0200aSDmitry Karpeev + n - number of values in the merged array (== an + bn) 611c1f0200aSDmitry Karpeev . I - merged sorted array 612c1f0200aSDmitry Karpeev - J - merged additional array 613c1f0200aSDmitry Karpeev 614c1f0200aSDmitry Karpeev Level: intermediate 615c1f0200aSDmitry Karpeev 616c1f0200aSDmitry Karpeev Concepts: merging^arrays 617c1f0200aSDmitry Karpeev 618c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 619c1f0200aSDmitry Karpeev @*/ 620d7aa01a8SHong Zhang PetscErrorCode PetscMergeIntArrayPair(PetscInt an,const PetscInt *aI, const PetscInt *aJ, PetscInt bn, const PetscInt *bI, const PetscInt *bJ, PetscInt *n, PetscInt **L, PetscInt **J) 621c1f0200aSDmitry Karpeev { 622c1f0200aSDmitry Karpeev PetscErrorCode ierr; 623d7aa01a8SHong Zhang PetscInt n_, *L_ = *L, *J_= *J, ak, bk, k; 624c1f0200aSDmitry Karpeev 625c1f0200aSDmitry Karpeev n_ = an + bn; 626c1f0200aSDmitry Karpeev *n = n_; 627d7aa01a8SHong Zhang if (!L_) { 628d7aa01a8SHong Zhang ierr = PetscMalloc(n_*sizeof(PetscInt), L); CHKERRQ(ierr); 629d7aa01a8SHong Zhang L_ = *L; 630c1f0200aSDmitry Karpeev } 631c1f0200aSDmitry Karpeev if (!J_){ 632c1f0200aSDmitry Karpeev ierr = PetscMalloc(n_*sizeof(PetscInt), &J_); CHKERRQ(ierr); 633c1f0200aSDmitry Karpeev J_ = *J; 634c1f0200aSDmitry Karpeev } 635c1f0200aSDmitry Karpeev k = ak = bk = 0; 636c1f0200aSDmitry Karpeev while(ak < an && bk < bn) { 637c1f0200aSDmitry Karpeev if (aI[ak] <= bI[bk]) { 638d7aa01a8SHong Zhang L_[k] = aI[ak]; 639c1f0200aSDmitry Karpeev J_[k] = aJ[ak]; 640c1f0200aSDmitry Karpeev ++ak; 641c1f0200aSDmitry Karpeev ++k; 642c1f0200aSDmitry Karpeev } 643c1f0200aSDmitry Karpeev else { 644d7aa01a8SHong Zhang L_[k] = bI[bk]; 645c1f0200aSDmitry Karpeev J_[k] = bJ[bk]; 646c1f0200aSDmitry Karpeev ++bk; 647c1f0200aSDmitry Karpeev ++k; 648c1f0200aSDmitry Karpeev } 649c1f0200aSDmitry Karpeev } 650c1f0200aSDmitry Karpeev if (ak < an) { 651d7aa01a8SHong Zhang ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt)); CHKERRQ(ierr); 652c1f0200aSDmitry Karpeev ierr = PetscMemcpy(J_+k,aJ+ak,(an-ak)*sizeof(PetscInt)); CHKERRQ(ierr); 653c1f0200aSDmitry Karpeev k += (an-ak); 654c1f0200aSDmitry Karpeev } 655c1f0200aSDmitry Karpeev if (bk < bn) { 656d7aa01a8SHong Zhang ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt)); CHKERRQ(ierr); 657c1f0200aSDmitry Karpeev ierr = PetscMemcpy(J_+k,bJ+bk,(bn-bk)*sizeof(PetscInt)); CHKERRQ(ierr); 658c1f0200aSDmitry Karpeev k += (bn-bk); 659c1f0200aSDmitry Karpeev } 660c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 661c1f0200aSDmitry Karpeev } 662c1f0200aSDmitry Karpeev 663e5c89e4eSSatish Balay 66442eaa7f3SBarry Smith #undef __FUNCT__ 66542eaa7f3SBarry Smith #define __FUNCT__ "PetscProcessTree" 66642eaa7f3SBarry Smith /*@ 66742eaa7f3SBarry Smith PetscProcessTree - Prepares tree data to be displayed graphically 66842eaa7f3SBarry Smith 66942eaa7f3SBarry Smith Not Collective 67042eaa7f3SBarry Smith 67142eaa7f3SBarry Smith Input Parameters: 67242eaa7f3SBarry Smith + n - number of values 67342eaa7f3SBarry Smith . mask - indicates those entries in the tree, location 0 is always masked 67442eaa7f3SBarry Smith - parentid - indicates the parent of each entry 67542eaa7f3SBarry Smith 67642eaa7f3SBarry Smith Output Parameters: 67742eaa7f3SBarry Smith + Nlevels - the number of levels 67842eaa7f3SBarry Smith . Level - for each node tells its level 67942eaa7f3SBarry Smith . Levelcnts - the number of nodes on each level 68042eaa7f3SBarry Smith . Idbylevel - a list of ids on each of the levels, first level followed by second etc 68142eaa7f3SBarry Smith - Column - for each id tells its column index 68242eaa7f3SBarry Smith 68342eaa7f3SBarry Smith Level: intermediate 68442eaa7f3SBarry Smith 68542eaa7f3SBarry Smith 68642eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation() 68742eaa7f3SBarry Smith @*/ 6887087cfbeSBarry Smith PetscErrorCode PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column) 68942eaa7f3SBarry Smith { 69042eaa7f3SBarry Smith PetscInt i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column; 69142eaa7f3SBarry Smith PetscErrorCode ierr; 692ace3abfcSBarry Smith PetscBool done = PETSC_FALSE; 69342eaa7f3SBarry Smith 69442eaa7f3SBarry Smith PetscFunctionBegin; 69542eaa7f3SBarry Smith if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set"); 69642eaa7f3SBarry Smith for (i=0; i<n; i++) { 69742eaa7f3SBarry Smith if (mask[i]) continue; 69842eaa7f3SBarry Smith if (parentid[i] == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent"); 69942eaa7f3SBarry Smith if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked"); 70042eaa7f3SBarry Smith } 70142eaa7f3SBarry Smith 70242eaa7f3SBarry Smith 70342eaa7f3SBarry Smith for (i=0; i<n; i++) { 70442eaa7f3SBarry Smith if (!mask[i]) nmask++; 70542eaa7f3SBarry Smith } 70642eaa7f3SBarry Smith 70742eaa7f3SBarry Smith /* determine the level in the tree of each node */ 70842eaa7f3SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&level);CHKERRQ(ierr); 70942eaa7f3SBarry Smith ierr = PetscMemzero(level,n*sizeof(PetscInt));CHKERRQ(ierr); 71042eaa7f3SBarry Smith level[0] = 1; 71142eaa7f3SBarry Smith while (!done) { 71242eaa7f3SBarry Smith done = PETSC_TRUE; 71342eaa7f3SBarry Smith for (i=0; i<n; i++) { 71442eaa7f3SBarry Smith if (mask[i]) continue; 71542eaa7f3SBarry Smith if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1; 71642eaa7f3SBarry Smith else if (!level[i]) done = PETSC_FALSE; 71742eaa7f3SBarry Smith } 71842eaa7f3SBarry Smith } 71942eaa7f3SBarry Smith for (i=0; i<n; i++) { 72042eaa7f3SBarry Smith level[i]--; 72142eaa7f3SBarry Smith nlevels = PetscMax(nlevels,level[i]); 72242eaa7f3SBarry Smith } 72342eaa7f3SBarry Smith 72442eaa7f3SBarry Smith /* count the number of nodes on each level and its max */ 72542eaa7f3SBarry Smith ierr = PetscMalloc(nlevels*sizeof(PetscInt),&levelcnt);CHKERRQ(ierr); 72642eaa7f3SBarry Smith ierr = PetscMemzero(levelcnt,nlevels*sizeof(PetscInt));CHKERRQ(ierr); 72742eaa7f3SBarry Smith for (i=0; i<n; i++) { 72842eaa7f3SBarry Smith if (mask[i]) continue; 72942eaa7f3SBarry Smith levelcnt[level[i]-1]++; 73042eaa7f3SBarry Smith } 73142eaa7f3SBarry Smith for (i=0; i<nlevels;i++) { 73242eaa7f3SBarry Smith levelmax = PetscMax(levelmax,levelcnt[i]); 73342eaa7f3SBarry Smith } 73442eaa7f3SBarry Smith 73542eaa7f3SBarry Smith /* for each level sort the ids by the parent id */ 73642eaa7f3SBarry Smith ierr = PetscMalloc2(levelmax,PetscInt,&workid,levelmax,PetscInt,&workparentid);CHKERRQ(ierr); 73742eaa7f3SBarry Smith ierr = PetscMalloc(nmask*sizeof(PetscInt),&idbylevel);CHKERRQ(ierr); 73842eaa7f3SBarry Smith for (j=1; j<=nlevels;j++) { 73942eaa7f3SBarry Smith cnt = 0; 74042eaa7f3SBarry Smith for (i=0; i<n; i++) { 74142eaa7f3SBarry Smith if (mask[i]) continue; 74242eaa7f3SBarry Smith if (level[i] != j) continue; 74342eaa7f3SBarry Smith workid[cnt] = i; 74442eaa7f3SBarry Smith workparentid[cnt++] = parentid[i]; 74542eaa7f3SBarry Smith } 74642eaa7f3SBarry Smith /* PetscIntView(cnt,workparentid,0); 74742eaa7f3SBarry Smith PetscIntView(cnt,workid,0); 74842eaa7f3SBarry Smith ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr); 74942eaa7f3SBarry Smith PetscIntView(cnt,workparentid,0); 75042eaa7f3SBarry Smith PetscIntView(cnt,workid,0);*/ 75142eaa7f3SBarry Smith ierr = PetscMemcpy(idbylevel+tcnt,workid,cnt*sizeof(PetscInt));CHKERRQ(ierr); 75242eaa7f3SBarry Smith tcnt += cnt; 75342eaa7f3SBarry Smith } 75442eaa7f3SBarry Smith if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes"); 75542eaa7f3SBarry Smith ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr); 75642eaa7f3SBarry Smith 75742eaa7f3SBarry Smith /* for each node list its column */ 75842eaa7f3SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&column);CHKERRQ(ierr); 75942eaa7f3SBarry Smith cnt = 0; 76042eaa7f3SBarry Smith for (j=0; j<nlevels; j++) { 76142eaa7f3SBarry Smith for (i=0; i<levelcnt[j]; i++) { 76242eaa7f3SBarry Smith column[idbylevel[cnt++]] = i; 76342eaa7f3SBarry Smith } 76442eaa7f3SBarry Smith } 76542eaa7f3SBarry Smith 76642eaa7f3SBarry Smith *Nlevels = nlevels; 76742eaa7f3SBarry Smith *Level = level; 76842eaa7f3SBarry Smith *Levelcnt = levelcnt; 76942eaa7f3SBarry Smith *Idbylevel = idbylevel; 77042eaa7f3SBarry Smith *Column = column; 77142eaa7f3SBarry Smith PetscFunctionReturn(0); 77242eaa7f3SBarry Smith } 773