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 1261db36a52SBarry Smith 127e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 128e5c89e4eSSatish Balay #define SWAP2(a,b,c,d,t) {t=a;a=b;b=t;t=c;c=d;d=t;} 129e5c89e4eSSatish Balay 130e5c89e4eSSatish Balay #undef __FUNCT__ 131e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithArray_Private" 132e5c89e4eSSatish Balay /* 133e5c89e4eSSatish Balay A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 134e5c89e4eSSatish Balay Assumes 0 origin for v, number of elements = right+1 (right is index of 135e5c89e4eSSatish Balay right-most member). 136e5c89e4eSSatish Balay */ 137e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithArray_Private(PetscInt *v,PetscInt *V,PetscInt right) 138e5c89e4eSSatish Balay { 139e5c89e4eSSatish Balay PetscErrorCode ierr; 140e5c89e4eSSatish Balay PetscInt i,vl,last,tmp; 141e5c89e4eSSatish Balay 142e5c89e4eSSatish Balay PetscFunctionBegin; 143e5c89e4eSSatish Balay if (right <= 1) { 144e5c89e4eSSatish Balay if (right == 1) { 145e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp); 146e5c89e4eSSatish Balay } 147e5c89e4eSSatish Balay PetscFunctionReturn(0); 148e5c89e4eSSatish Balay } 149e5c89e4eSSatish Balay SWAP2(v[0],v[right/2],V[0],V[right/2],tmp); 150e5c89e4eSSatish Balay vl = v[0]; 151e5c89e4eSSatish Balay last = 0; 152e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 153e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);} 154e5c89e4eSSatish Balay } 155e5c89e4eSSatish Balay SWAP2(v[0],v[last],V[0],V[last],tmp); 156e5c89e4eSSatish Balay ierr = PetscSortIntWithArray_Private(v,V,last-1);CHKERRQ(ierr); 157e5c89e4eSSatish Balay ierr = PetscSortIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 158e5c89e4eSSatish Balay PetscFunctionReturn(0); 159e5c89e4eSSatish Balay } 160e5c89e4eSSatish Balay 161e5c89e4eSSatish Balay #undef __FUNCT__ 162e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithArray" 163e5c89e4eSSatish Balay /*@ 164e5c89e4eSSatish Balay PetscSortIntWithArray - Sorts an array of integers in place in increasing order; 165e5c89e4eSSatish Balay changes a second array to match the sorted first array. 166e5c89e4eSSatish Balay 167e5c89e4eSSatish Balay Not Collective 168e5c89e4eSSatish Balay 169e5c89e4eSSatish Balay Input Parameters: 170e5c89e4eSSatish Balay + n - number of values 171e5c89e4eSSatish Balay . i - array of integers 172e5c89e4eSSatish Balay - I - second array of integers 173e5c89e4eSSatish Balay 174e5c89e4eSSatish Balay Level: intermediate 175e5c89e4eSSatish Balay 176e5c89e4eSSatish Balay Concepts: sorting^ints with array 177e5c89e4eSSatish Balay 178e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 179e5c89e4eSSatish Balay @*/ 1807087cfbeSBarry Smith PetscErrorCode PetscSortIntWithArray(PetscInt n,PetscInt i[],PetscInt Ii[]) 181e5c89e4eSSatish Balay { 182e5c89e4eSSatish Balay PetscErrorCode ierr; 183e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 184e5c89e4eSSatish Balay 185e5c89e4eSSatish Balay PetscFunctionBegin; 186e5c89e4eSSatish Balay if (n<8) { 187e5c89e4eSSatish Balay for (k=0; k<n; k++) { 188e5c89e4eSSatish Balay ik = i[k]; 189e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 190e5c89e4eSSatish Balay if (ik > i[j]) { 191b7940d39SSatish Balay SWAP2(i[k],i[j],Ii[k],Ii[j],tmp); 192e5c89e4eSSatish Balay ik = i[k]; 193e5c89e4eSSatish Balay } 194e5c89e4eSSatish Balay } 195e5c89e4eSSatish Balay } 196e5c89e4eSSatish Balay } else { 197b7940d39SSatish Balay ierr = PetscSortIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr); 198e5c89e4eSSatish Balay } 199e5c89e4eSSatish Balay PetscFunctionReturn(0); 200e5c89e4eSSatish Balay } 201e5c89e4eSSatish Balay 202c1f0200aSDmitry Karpeev 203c1f0200aSDmitry 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;} 204c1f0200aSDmitry Karpeev 205c1f0200aSDmitry Karpeev #undef __FUNCT__ 206c1f0200aSDmitry Karpeev #define __FUNCT__ "PetscSortIntWithArrayPair_Private" 207c1f0200aSDmitry Karpeev /* 208c1f0200aSDmitry Karpeev A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 209c1f0200aSDmitry Karpeev Assumes 0 origin for v, number of elements = right+1 (right is index of 210c1f0200aSDmitry Karpeev right-most member). 211c1f0200aSDmitry Karpeev */ 212d7aa01a8SHong Zhang static PetscErrorCode PetscSortIntWithArrayPair_Private(PetscInt *L,PetscInt *J, PetscInt *K, PetscInt right) 213c1f0200aSDmitry Karpeev { 214c1f0200aSDmitry Karpeev PetscErrorCode ierr; 215c1f0200aSDmitry Karpeev PetscInt i,vl,last,tmp; 216c1f0200aSDmitry Karpeev 217c1f0200aSDmitry Karpeev PetscFunctionBegin; 218c1f0200aSDmitry Karpeev if (right <= 1) { 219c1f0200aSDmitry Karpeev if (right == 1) { 220d7aa01a8SHong Zhang if (L[0] > L[1]) SWAP3(L[0],L[1],J[0],J[1],K[0],K[1],tmp); 221c1f0200aSDmitry Karpeev } 222c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 223c1f0200aSDmitry Karpeev } 224d7aa01a8SHong Zhang SWAP3(L[0],L[right/2],J[0],J[right/2],K[0],K[right/2],tmp); 225d7aa01a8SHong Zhang vl = L[0]; 226c1f0200aSDmitry Karpeev last = 0; 227c1f0200aSDmitry Karpeev for (i=1; i<=right; i++) { 228d7aa01a8SHong Zhang if (L[i] < vl) {last++; SWAP3(L[last],L[i],J[last],J[i],K[last],K[i],tmp);} 229c1f0200aSDmitry Karpeev } 230d7aa01a8SHong Zhang SWAP3(L[0],L[last],J[0],J[last],K[0],K[last],tmp); 231d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L,J,K,last-1);CHKERRQ(ierr); 232d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L+last+1,J+last+1,K+last+1,right-(last+1));CHKERRQ(ierr); 233c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 234c1f0200aSDmitry Karpeev } 235c1f0200aSDmitry Karpeev 236c1f0200aSDmitry Karpeev #undef __FUNCT__ 237c1f0200aSDmitry Karpeev #define __FUNCT__ "PetscSortIntWithArrayPair" 238c1f0200aSDmitry Karpeev /*@ 239c1f0200aSDmitry Karpeev PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order; 240c1f0200aSDmitry Karpeev changes a pair of integer arrays to match the sorted first array. 241c1f0200aSDmitry Karpeev 242c1f0200aSDmitry Karpeev Not Collective 243c1f0200aSDmitry Karpeev 244c1f0200aSDmitry Karpeev Input Parameters: 245c1f0200aSDmitry Karpeev + n - number of values 246c1f0200aSDmitry Karpeev . I - array of integers 247c1f0200aSDmitry Karpeev . J - second array of integers (first array of the pair) 248c1f0200aSDmitry Karpeev - K - third array of integers (second array of the pair) 249c1f0200aSDmitry Karpeev 250c1f0200aSDmitry Karpeev Level: intermediate 251c1f0200aSDmitry Karpeev 252c1f0200aSDmitry Karpeev Concepts: sorting^ints with array pair 253c1f0200aSDmitry Karpeev 254c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray() 255c1f0200aSDmitry Karpeev @*/ 256d7aa01a8SHong Zhang PetscErrorCode PetscSortIntWithArrayPair(PetscInt n,PetscInt *L,PetscInt *J, PetscInt *K) 257c1f0200aSDmitry Karpeev { 258c1f0200aSDmitry Karpeev PetscErrorCode ierr; 259c1f0200aSDmitry Karpeev PetscInt j,k,tmp,ik; 260c1f0200aSDmitry Karpeev 261c1f0200aSDmitry Karpeev PetscFunctionBegin; 262c1f0200aSDmitry Karpeev if (n<8) { 263c1f0200aSDmitry Karpeev for (k=0; k<n; k++) { 264d7aa01a8SHong Zhang ik = L[k]; 265c1f0200aSDmitry Karpeev for (j=k+1; j<n; j++) { 266d7aa01a8SHong Zhang if (ik > L[j]) { 267d7aa01a8SHong Zhang SWAP3(L[k],L[j],J[k],J[j],K[k],K[j],tmp); 268d7aa01a8SHong Zhang ik = L[k]; 269c1f0200aSDmitry Karpeev } 270c1f0200aSDmitry Karpeev } 271c1f0200aSDmitry Karpeev } 272c1f0200aSDmitry Karpeev } else { 273d7aa01a8SHong Zhang ierr = PetscSortIntWithArrayPair_Private(L,J,K,n-1);CHKERRQ(ierr); 274c1f0200aSDmitry Karpeev } 275c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 276c1f0200aSDmitry Karpeev } 277c1f0200aSDmitry Karpeev 278*17d7d925SStefano Zampini #undef __FUNCT__ 279*17d7d925SStefano Zampini #define __FUNCT__ "PetscSortMPIInt_Private" 280*17d7d925SStefano Zampini /* 281*17d7d925SStefano Zampini A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 282*17d7d925SStefano Zampini Assumes 0 origin for v, number of elements = right+1 (right is index of 283*17d7d925SStefano Zampini right-most member). 284*17d7d925SStefano Zampini */ 285*17d7d925SStefano Zampini static void PetscSortMPIInt_Private(PetscMPIInt *v,PetscInt right) 286*17d7d925SStefano Zampini { 287*17d7d925SStefano Zampini PetscInt i,j; 288*17d7d925SStefano Zampini PetscMPIInt pivot,tmp; 289*17d7d925SStefano Zampini 290*17d7d925SStefano Zampini if (right <= 1) { 291*17d7d925SStefano Zampini if (right == 1) { 292*17d7d925SStefano Zampini if (v[0] > v[1]) SWAP(v[0],v[1],tmp); 293*17d7d925SStefano Zampini } 294*17d7d925SStefano Zampini return; 295*17d7d925SStefano Zampini } 296*17d7d925SStefano Zampini i = MEDIAN(v,right); /* Choose a pivot */ 297*17d7d925SStefano Zampini SWAP(v[0],v[i],tmp); /* Move it out of the way */ 298*17d7d925SStefano Zampini pivot = v[0]; 299*17d7d925SStefano Zampini for (i=0,j=right+1;;) { 300*17d7d925SStefano Zampini while (++i < j && v[i] <= pivot); /* Scan from the left */ 301*17d7d925SStefano Zampini while (v[--j] > pivot); /* Scan from the right */ 302*17d7d925SStefano Zampini if (i >= j) break; 303*17d7d925SStefano Zampini SWAP(v[i],v[j],tmp); 304*17d7d925SStefano Zampini } 305*17d7d925SStefano Zampini SWAP(v[0],v[j],tmp); /* Put pivot back in place. */ 306*17d7d925SStefano Zampini PetscSortMPIInt_Private(v,j-1); 307*17d7d925SStefano Zampini PetscSortMPIInt_Private(v+j+1,right-(j+1)); 308*17d7d925SStefano Zampini } 309*17d7d925SStefano Zampini 310*17d7d925SStefano Zampini #undef __FUNCT__ 311*17d7d925SStefano Zampini #define __FUNCT__ "PetscSortMPIInt" 312*17d7d925SStefano Zampini /*@ 313*17d7d925SStefano Zampini PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order. 314*17d7d925SStefano Zampini 315*17d7d925SStefano Zampini Not Collective 316*17d7d925SStefano Zampini 317*17d7d925SStefano Zampini Input Parameters: 318*17d7d925SStefano Zampini + n - number of values 319*17d7d925SStefano Zampini - i - array of integers 320*17d7d925SStefano Zampini 321*17d7d925SStefano Zampini Level: intermediate 322*17d7d925SStefano Zampini 323*17d7d925SStefano Zampini Concepts: sorting^ints 324*17d7d925SStefano Zampini 325*17d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation() 326*17d7d925SStefano Zampini @*/ 327*17d7d925SStefano Zampini PetscErrorCode PetscSortMPIInt(PetscInt n,PetscMPIInt i[]) 328*17d7d925SStefano Zampini { 329*17d7d925SStefano Zampini PetscInt j,k; 330*17d7d925SStefano Zampini PetscMPIInt tmp,ik; 331*17d7d925SStefano Zampini 332*17d7d925SStefano Zampini PetscFunctionBegin; 333*17d7d925SStefano Zampini if (n<8) { 334*17d7d925SStefano Zampini for (k=0; k<n; k++) { 335*17d7d925SStefano Zampini ik = i[k]; 336*17d7d925SStefano Zampini for (j=k+1; j<n; j++) { 337*17d7d925SStefano Zampini if (ik > i[j]) { 338*17d7d925SStefano Zampini SWAP(i[k],i[j],tmp); 339*17d7d925SStefano Zampini ik = i[k]; 340*17d7d925SStefano Zampini } 341*17d7d925SStefano Zampini } 342*17d7d925SStefano Zampini } 343*17d7d925SStefano Zampini } else { 344*17d7d925SStefano Zampini PetscSortMPIInt_Private(i,n-1); 345*17d7d925SStefano Zampini } 346*17d7d925SStefano Zampini PetscFunctionReturn(0); 347*17d7d925SStefano Zampini } 348*17d7d925SStefano Zampini 349*17d7d925SStefano Zampini #undef __FUNCT__ 350*17d7d925SStefano Zampini #define __FUNCT__ "PetscSortRemoveDupsMPIInt" 351*17d7d925SStefano Zampini /*@ 352*17d7d925SStefano Zampini PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries 353*17d7d925SStefano Zampini 354*17d7d925SStefano Zampini Not Collective 355*17d7d925SStefano Zampini 356*17d7d925SStefano Zampini Input Parameters: 357*17d7d925SStefano Zampini + n - number of values 358*17d7d925SStefano Zampini - ii - array of integers 359*17d7d925SStefano Zampini 360*17d7d925SStefano Zampini Output Parameter: 361*17d7d925SStefano Zampini . n - number of non-redundant values 362*17d7d925SStefano Zampini 363*17d7d925SStefano Zampini Level: intermediate 364*17d7d925SStefano Zampini 365*17d7d925SStefano Zampini Concepts: sorting^ints 366*17d7d925SStefano Zampini 367*17d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt() 368*17d7d925SStefano Zampini @*/ 369*17d7d925SStefano Zampini PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt ii[]) 370*17d7d925SStefano Zampini { 371*17d7d925SStefano Zampini PetscErrorCode ierr; 372*17d7d925SStefano Zampini PetscInt i,s = 0,N = *n, b = 0; 373*17d7d925SStefano Zampini 374*17d7d925SStefano Zampini PetscFunctionBegin; 375*17d7d925SStefano Zampini ierr = PetscSortMPIInt(N,ii);CHKERRQ(ierr); 376*17d7d925SStefano Zampini for (i=0; i<N-1; i++) { 377*17d7d925SStefano Zampini if (ii[b+s+1] != ii[b]) {ii[b+1] = ii[b+s+1]; b++;} 378*17d7d925SStefano Zampini else s++; 379*17d7d925SStefano Zampini } 380*17d7d925SStefano Zampini *n = N - s; 381*17d7d925SStefano Zampini PetscFunctionReturn(0); 382*17d7d925SStefano Zampini } 383c1f0200aSDmitry Karpeev 3844d615ea0SBarry Smith #undef __FUNCT__ 3854d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray_Private" 3864d615ea0SBarry Smith /* 3874d615ea0SBarry Smith A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 3884d615ea0SBarry Smith Assumes 0 origin for v, number of elements = right+1 (right is index of 3894d615ea0SBarry Smith right-most member). 3904d615ea0SBarry Smith */ 3914d615ea0SBarry Smith static PetscErrorCode PetscSortMPIIntWithArray_Private(PetscMPIInt *v,PetscMPIInt *V,PetscMPIInt right) 3924d615ea0SBarry Smith { 3934d615ea0SBarry Smith PetscErrorCode ierr; 3944d615ea0SBarry Smith PetscMPIInt i,vl,last,tmp; 3954d615ea0SBarry Smith 3964d615ea0SBarry Smith PetscFunctionBegin; 3974d615ea0SBarry Smith if (right <= 1) { 3984d615ea0SBarry Smith if (right == 1) { 3994d615ea0SBarry Smith if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp); 4004d615ea0SBarry Smith } 4014d615ea0SBarry Smith PetscFunctionReturn(0); 4024d615ea0SBarry Smith } 4034d615ea0SBarry Smith SWAP2(v[0],v[right/2],V[0],V[right/2],tmp); 4044d615ea0SBarry Smith vl = v[0]; 4054d615ea0SBarry Smith last = 0; 4064d615ea0SBarry Smith for (i=1; i<=right; i++) { 4074d615ea0SBarry Smith if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);} 4084d615ea0SBarry Smith } 4094d615ea0SBarry Smith SWAP2(v[0],v[last],V[0],V[last],tmp); 4104d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v,V,last-1);CHKERRQ(ierr); 4114d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 4124d615ea0SBarry Smith PetscFunctionReturn(0); 4134d615ea0SBarry Smith } 4144d615ea0SBarry Smith 4154d615ea0SBarry Smith #undef __FUNCT__ 4164d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray" 4174d615ea0SBarry Smith /*@ 4184d615ea0SBarry Smith PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order; 4194d615ea0SBarry Smith changes a second array to match the sorted first array. 4204d615ea0SBarry Smith 4214d615ea0SBarry Smith Not Collective 4224d615ea0SBarry Smith 4234d615ea0SBarry Smith Input Parameters: 4244d615ea0SBarry Smith + n - number of values 4254d615ea0SBarry Smith . i - array of integers 4264d615ea0SBarry Smith - I - second array of integers 4274d615ea0SBarry Smith 4284d615ea0SBarry Smith Level: intermediate 4294d615ea0SBarry Smith 4304d615ea0SBarry Smith Concepts: sorting^ints with array 4314d615ea0SBarry Smith 4324d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 4334d615ea0SBarry Smith @*/ 4347087cfbeSBarry Smith PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt i[],PetscMPIInt Ii[]) 4354d615ea0SBarry Smith { 4364d615ea0SBarry Smith PetscErrorCode ierr; 4374d615ea0SBarry Smith PetscMPIInt j,k,tmp,ik; 4384d615ea0SBarry Smith 4394d615ea0SBarry Smith PetscFunctionBegin; 4404d615ea0SBarry Smith if (n<8) { 4414d615ea0SBarry Smith for (k=0; k<n; k++) { 4424d615ea0SBarry Smith ik = i[k]; 4434d615ea0SBarry Smith for (j=k+1; j<n; j++) { 4444d615ea0SBarry Smith if (ik > i[j]) { 4454d615ea0SBarry Smith SWAP2(i[k],i[j],Ii[k],Ii[j],tmp); 4464d615ea0SBarry Smith ik = i[k]; 4474d615ea0SBarry Smith } 4484d615ea0SBarry Smith } 4494d615ea0SBarry Smith } 4504d615ea0SBarry Smith } else { 4514d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr); 4524d615ea0SBarry Smith } 4534d615ea0SBarry Smith PetscFunctionReturn(0); 4544d615ea0SBarry Smith } 4554d615ea0SBarry Smith 456e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 457e5c89e4eSSatish Balay #define SWAP2IntScalar(a,b,c,d,t,ts) {t=a;a=b;b=t;ts=c;c=d;d=ts;} 458e5c89e4eSSatish Balay 459e5c89e4eSSatish Balay #undef __FUNCT__ 460e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray_Private" 461e5c89e4eSSatish Balay /* 462e5c89e4eSSatish Balay Modified from PetscSortIntWithArray_Private(). 463e5c89e4eSSatish Balay */ 464e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithScalarArray_Private(PetscInt *v,PetscScalar *V,PetscInt right) 465e5c89e4eSSatish Balay { 466e5c89e4eSSatish Balay PetscErrorCode ierr; 467e5c89e4eSSatish Balay PetscInt i,vl,last,tmp; 468e5c89e4eSSatish Balay PetscScalar stmp; 469e5c89e4eSSatish Balay 470e5c89e4eSSatish Balay PetscFunctionBegin; 471e5c89e4eSSatish Balay if (right <= 1) { 472e5c89e4eSSatish Balay if (right == 1) { 473e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP2IntScalar(v[0],v[1],V[0],V[1],tmp,stmp); 474e5c89e4eSSatish Balay } 475e5c89e4eSSatish Balay PetscFunctionReturn(0); 476e5c89e4eSSatish Balay } 477e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[right/2],V[0],V[right/2],tmp,stmp); 478e5c89e4eSSatish Balay vl = v[0]; 479e5c89e4eSSatish Balay last = 0; 480e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 481e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP2IntScalar(v[last],v[i],V[last],V[i],tmp,stmp);} 482e5c89e4eSSatish Balay } 483e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[last],V[0],V[last],tmp,stmp); 484e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v,V,last-1);CHKERRQ(ierr); 485e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 486e5c89e4eSSatish Balay PetscFunctionReturn(0); 487e5c89e4eSSatish Balay } 488e5c89e4eSSatish Balay 489e5c89e4eSSatish Balay #undef __FUNCT__ 490e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray" 491e5c89e4eSSatish Balay /*@ 492e5c89e4eSSatish Balay PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order; 493e5c89e4eSSatish Balay changes a second SCALAR array to match the sorted first INTEGER array. 494e5c89e4eSSatish Balay 495e5c89e4eSSatish Balay Not Collective 496e5c89e4eSSatish Balay 497e5c89e4eSSatish Balay Input Parameters: 498e5c89e4eSSatish Balay + n - number of values 499e5c89e4eSSatish Balay . i - array of integers 500e5c89e4eSSatish Balay - I - second array of scalars 501e5c89e4eSSatish Balay 502e5c89e4eSSatish Balay Level: intermediate 503e5c89e4eSSatish Balay 504e5c89e4eSSatish Balay Concepts: sorting^ints with array 505e5c89e4eSSatish Balay 506e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 507e5c89e4eSSatish Balay @*/ 5087087cfbeSBarry Smith PetscErrorCode PetscSortIntWithScalarArray(PetscInt n,PetscInt i[],PetscScalar Ii[]) 509e5c89e4eSSatish Balay { 510e5c89e4eSSatish Balay PetscErrorCode ierr; 511e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 512e5c89e4eSSatish Balay PetscScalar stmp; 513e5c89e4eSSatish Balay 514e5c89e4eSSatish Balay PetscFunctionBegin; 515e5c89e4eSSatish Balay if (n<8) { 516e5c89e4eSSatish Balay for (k=0; k<n; k++) { 517e5c89e4eSSatish Balay ik = i[k]; 518e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 519e5c89e4eSSatish Balay if (ik > i[j]) { 520b7940d39SSatish Balay SWAP2IntScalar(i[k],i[j],Ii[k],Ii[j],tmp,stmp); 521e5c89e4eSSatish Balay ik = i[k]; 522e5c89e4eSSatish Balay } 523e5c89e4eSSatish Balay } 524e5c89e4eSSatish Balay } 525e5c89e4eSSatish Balay } else { 526b7940d39SSatish Balay ierr = PetscSortIntWithScalarArray_Private(i,Ii,n-1);CHKERRQ(ierr); 527e5c89e4eSSatish Balay } 528e5c89e4eSSatish Balay PetscFunctionReturn(0); 529e5c89e4eSSatish Balay } 530e5c89e4eSSatish Balay 531b4301105SBarry Smith 532b4301105SBarry Smith 533b4301105SBarry Smith #undef __FUNCT__ 534c1f0200aSDmitry Karpeev #define __FUNCT__ "PetscMergeIntArrayPair" 535c1f0200aSDmitry Karpeev /*@ 536c1f0200aSDmitry Karpeev PetscMergeIntArrayPair - Merges two SORTED integer arrays along with an additional array of integers. 537c1f0200aSDmitry Karpeev The additional arrays are the same length as sorted arrays and are merged 538c1f0200aSDmitry Karpeev in the order determined by the merging of the sorted pair. 539c1f0200aSDmitry Karpeev 540c1f0200aSDmitry Karpeev Not Collective 541c1f0200aSDmitry Karpeev 542c1f0200aSDmitry Karpeev Input Parameters: 543c1f0200aSDmitry Karpeev + an - number of values in the first array 544c1f0200aSDmitry Karpeev . aI - first sorted array of integers 545c1f0200aSDmitry Karpeev . aJ - first additional array of integers 546c1f0200aSDmitry Karpeev . bn - number of values in the second array 547c1f0200aSDmitry Karpeev . bI - second array of integers 548c1f0200aSDmitry Karpeev - bJ - second additional of integers 549c1f0200aSDmitry Karpeev 550c1f0200aSDmitry Karpeev Output Parameters: 551c1f0200aSDmitry Karpeev + n - number of values in the merged array (== an + bn) 552c1f0200aSDmitry Karpeev . I - merged sorted array 553c1f0200aSDmitry Karpeev - J - merged additional array 554c1f0200aSDmitry Karpeev 555c1f0200aSDmitry Karpeev Level: intermediate 556c1f0200aSDmitry Karpeev 557c1f0200aSDmitry Karpeev Concepts: merging^arrays 558c1f0200aSDmitry Karpeev 559c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 560c1f0200aSDmitry Karpeev @*/ 561d7aa01a8SHong 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) 562c1f0200aSDmitry Karpeev { 563c1f0200aSDmitry Karpeev PetscErrorCode ierr; 564d7aa01a8SHong Zhang PetscInt n_, *L_ = *L, *J_= *J, ak, bk, k; 565c1f0200aSDmitry Karpeev 566c1f0200aSDmitry Karpeev n_ = an + bn; 567c1f0200aSDmitry Karpeev *n = n_; 568d7aa01a8SHong Zhang if(!L_) { 569d7aa01a8SHong Zhang ierr = PetscMalloc(n_*sizeof(PetscInt), L); CHKERRQ(ierr); 570d7aa01a8SHong Zhang L_ = *L; 571c1f0200aSDmitry Karpeev } 572c1f0200aSDmitry Karpeev if(!J_){ 573c1f0200aSDmitry Karpeev ierr = PetscMalloc(n_*sizeof(PetscInt), &J_); CHKERRQ(ierr); 574c1f0200aSDmitry Karpeev J_ = *J; 575c1f0200aSDmitry Karpeev } 576c1f0200aSDmitry Karpeev k = ak = bk = 0; 577c1f0200aSDmitry Karpeev while(ak < an && bk < bn) { 578c1f0200aSDmitry Karpeev if(aI[ak] <= bI[bk]) { 579d7aa01a8SHong Zhang L_[k] = aI[ak]; 580c1f0200aSDmitry Karpeev J_[k] = aJ[ak]; 581c1f0200aSDmitry Karpeev ++ak; 582c1f0200aSDmitry Karpeev ++k; 583c1f0200aSDmitry Karpeev } 584c1f0200aSDmitry Karpeev else { 585d7aa01a8SHong Zhang L_[k] = bI[bk]; 586c1f0200aSDmitry Karpeev J_[k] = bJ[bk]; 587c1f0200aSDmitry Karpeev ++bk; 588c1f0200aSDmitry Karpeev ++k; 589c1f0200aSDmitry Karpeev } 590c1f0200aSDmitry Karpeev } 591c1f0200aSDmitry Karpeev if(ak < an) { 592d7aa01a8SHong Zhang ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt)); CHKERRQ(ierr); 593c1f0200aSDmitry Karpeev ierr = PetscMemcpy(J_+k,aJ+ak,(an-ak)*sizeof(PetscInt)); CHKERRQ(ierr); 594c1f0200aSDmitry Karpeev k += (an-ak); 595c1f0200aSDmitry Karpeev } 596c1f0200aSDmitry Karpeev if(bk < bn) { 597d7aa01a8SHong Zhang ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt)); CHKERRQ(ierr); 598c1f0200aSDmitry Karpeev ierr = PetscMemcpy(J_+k,bJ+bk,(bn-bk)*sizeof(PetscInt)); CHKERRQ(ierr); 599c1f0200aSDmitry Karpeev k += (bn-bk); 600c1f0200aSDmitry Karpeev } 601c1f0200aSDmitry Karpeev PetscFunctionReturn(0); 602c1f0200aSDmitry Karpeev } 603c1f0200aSDmitry Karpeev 604e5c89e4eSSatish Balay 60542eaa7f3SBarry Smith #undef __FUNCT__ 60642eaa7f3SBarry Smith #define __FUNCT__ "PetscProcessTree" 60742eaa7f3SBarry Smith /*@ 60842eaa7f3SBarry Smith PetscProcessTree - Prepares tree data to be displayed graphically 60942eaa7f3SBarry Smith 61042eaa7f3SBarry Smith Not Collective 61142eaa7f3SBarry Smith 61242eaa7f3SBarry Smith Input Parameters: 61342eaa7f3SBarry Smith + n - number of values 61442eaa7f3SBarry Smith . mask - indicates those entries in the tree, location 0 is always masked 61542eaa7f3SBarry Smith - parentid - indicates the parent of each entry 61642eaa7f3SBarry Smith 61742eaa7f3SBarry Smith Output Parameters: 61842eaa7f3SBarry Smith + Nlevels - the number of levels 61942eaa7f3SBarry Smith . Level - for each node tells its level 62042eaa7f3SBarry Smith . Levelcnts - the number of nodes on each level 62142eaa7f3SBarry Smith . Idbylevel - a list of ids on each of the levels, first level followed by second etc 62242eaa7f3SBarry Smith - Column - for each id tells its column index 62342eaa7f3SBarry Smith 62442eaa7f3SBarry Smith Level: intermediate 62542eaa7f3SBarry Smith 62642eaa7f3SBarry Smith 62742eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation() 62842eaa7f3SBarry Smith @*/ 6297087cfbeSBarry Smith PetscErrorCode PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column) 63042eaa7f3SBarry Smith { 63142eaa7f3SBarry Smith PetscInt i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column; 63242eaa7f3SBarry Smith PetscErrorCode ierr; 633ace3abfcSBarry Smith PetscBool done = PETSC_FALSE; 63442eaa7f3SBarry Smith 63542eaa7f3SBarry Smith PetscFunctionBegin; 63642eaa7f3SBarry Smith if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set"); 63742eaa7f3SBarry Smith for (i=0; i<n; i++) { 63842eaa7f3SBarry Smith if (mask[i]) continue; 63942eaa7f3SBarry Smith if (parentid[i] == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent"); 64042eaa7f3SBarry Smith if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked"); 64142eaa7f3SBarry Smith } 64242eaa7f3SBarry Smith 64342eaa7f3SBarry Smith 64442eaa7f3SBarry Smith for (i=0; i<n; i++) { 64542eaa7f3SBarry Smith if (!mask[i]) nmask++; 64642eaa7f3SBarry Smith } 64742eaa7f3SBarry Smith 64842eaa7f3SBarry Smith /* determine the level in the tree of each node */ 64942eaa7f3SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&level);CHKERRQ(ierr); 65042eaa7f3SBarry Smith ierr = PetscMemzero(level,n*sizeof(PetscInt));CHKERRQ(ierr); 65142eaa7f3SBarry Smith level[0] = 1; 65242eaa7f3SBarry Smith while (!done) { 65342eaa7f3SBarry Smith done = PETSC_TRUE; 65442eaa7f3SBarry Smith for (i=0; i<n; i++) { 65542eaa7f3SBarry Smith if (mask[i]) continue; 65642eaa7f3SBarry Smith if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1; 65742eaa7f3SBarry Smith else if (!level[i]) done = PETSC_FALSE; 65842eaa7f3SBarry Smith } 65942eaa7f3SBarry Smith } 66042eaa7f3SBarry Smith for (i=0; i<n; i++) { 66142eaa7f3SBarry Smith level[i]--; 66242eaa7f3SBarry Smith nlevels = PetscMax(nlevels,level[i]); 66342eaa7f3SBarry Smith } 66442eaa7f3SBarry Smith 66542eaa7f3SBarry Smith /* count the number of nodes on each level and its max */ 66642eaa7f3SBarry Smith ierr = PetscMalloc(nlevels*sizeof(PetscInt),&levelcnt);CHKERRQ(ierr); 66742eaa7f3SBarry Smith ierr = PetscMemzero(levelcnt,nlevels*sizeof(PetscInt));CHKERRQ(ierr); 66842eaa7f3SBarry Smith for (i=0; i<n; i++) { 66942eaa7f3SBarry Smith if (mask[i]) continue; 67042eaa7f3SBarry Smith levelcnt[level[i]-1]++; 67142eaa7f3SBarry Smith } 67242eaa7f3SBarry Smith for (i=0; i<nlevels;i++) { 67342eaa7f3SBarry Smith levelmax = PetscMax(levelmax,levelcnt[i]); 67442eaa7f3SBarry Smith } 67542eaa7f3SBarry Smith 67642eaa7f3SBarry Smith /* for each level sort the ids by the parent id */ 67742eaa7f3SBarry Smith ierr = PetscMalloc2(levelmax,PetscInt,&workid,levelmax,PetscInt,&workparentid);CHKERRQ(ierr); 67842eaa7f3SBarry Smith ierr = PetscMalloc(nmask*sizeof(PetscInt),&idbylevel);CHKERRQ(ierr); 67942eaa7f3SBarry Smith for (j=1; j<=nlevels;j++) { 68042eaa7f3SBarry Smith cnt = 0; 68142eaa7f3SBarry Smith for (i=0; i<n; i++) { 68242eaa7f3SBarry Smith if (mask[i]) continue; 68342eaa7f3SBarry Smith if (level[i] != j) continue; 68442eaa7f3SBarry Smith workid[cnt] = i; 68542eaa7f3SBarry Smith workparentid[cnt++] = parentid[i]; 68642eaa7f3SBarry Smith } 68742eaa7f3SBarry Smith /* PetscIntView(cnt,workparentid,0); 68842eaa7f3SBarry Smith PetscIntView(cnt,workid,0); 68942eaa7f3SBarry Smith ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr); 69042eaa7f3SBarry Smith PetscIntView(cnt,workparentid,0); 69142eaa7f3SBarry Smith PetscIntView(cnt,workid,0);*/ 69242eaa7f3SBarry Smith ierr = PetscMemcpy(idbylevel+tcnt,workid,cnt*sizeof(PetscInt));CHKERRQ(ierr); 69342eaa7f3SBarry Smith tcnt += cnt; 69442eaa7f3SBarry Smith } 69542eaa7f3SBarry Smith if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes"); 69642eaa7f3SBarry Smith ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr); 69742eaa7f3SBarry Smith 69842eaa7f3SBarry Smith /* for each node list its column */ 69942eaa7f3SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&column);CHKERRQ(ierr); 70042eaa7f3SBarry Smith cnt = 0; 70142eaa7f3SBarry Smith for (j=0; j<nlevels; j++) { 70242eaa7f3SBarry Smith for (i=0; i<levelcnt[j]; i++) { 70342eaa7f3SBarry Smith column[idbylevel[cnt++]] = i; 70442eaa7f3SBarry Smith } 70542eaa7f3SBarry Smith } 70642eaa7f3SBarry Smith 70742eaa7f3SBarry Smith *Nlevels = nlevels; 70842eaa7f3SBarry Smith *Level = level; 70942eaa7f3SBarry Smith *Levelcnt = levelcnt; 71042eaa7f3SBarry Smith *Idbylevel = idbylevel; 71142eaa7f3SBarry Smith *Column = column; 71242eaa7f3SBarry Smith PetscFunctionReturn(0); 71342eaa7f3SBarry Smith } 714