1e5c89e4eSSatish Balay #define PETSC_DLL 2e5c89e4eSSatish Balay /* 3e5c89e4eSSatish Balay This file contains routines for sorting integers. Values are sorted in place. 4e5c89e4eSSatish Balay */ 5d382aafbSBarry Smith #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 @*/ 708738c821SJed Brown PetscErrorCode PETSCSYS_DLLEXPORT 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 @*/ 1118738c821SJed Brown PetscErrorCode PETSCSYS_DLLEXPORT 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 @*/ 1808738c821SJed Brown PetscErrorCode PETSCSYS_DLLEXPORT 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 2024d615ea0SBarry Smith #undef __FUNCT__ 2034d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray_Private" 2044d615ea0SBarry Smith /* 2054d615ea0SBarry Smith A simple version of quicksort; taken from Kernighan and Ritchie, page 87. 2064d615ea0SBarry Smith Assumes 0 origin for v, number of elements = right+1 (right is index of 2074d615ea0SBarry Smith right-most member). 2084d615ea0SBarry Smith */ 2094d615ea0SBarry Smith static PetscErrorCode PetscSortMPIIntWithArray_Private(PetscMPIInt *v,PetscMPIInt *V,PetscMPIInt right) 2104d615ea0SBarry Smith { 2114d615ea0SBarry Smith PetscErrorCode ierr; 2124d615ea0SBarry Smith PetscMPIInt i,vl,last,tmp; 2134d615ea0SBarry Smith 2144d615ea0SBarry Smith PetscFunctionBegin; 2154d615ea0SBarry Smith if (right <= 1) { 2164d615ea0SBarry Smith if (right == 1) { 2174d615ea0SBarry Smith if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp); 2184d615ea0SBarry Smith } 2194d615ea0SBarry Smith PetscFunctionReturn(0); 2204d615ea0SBarry Smith } 2214d615ea0SBarry Smith SWAP2(v[0],v[right/2],V[0],V[right/2],tmp); 2224d615ea0SBarry Smith vl = v[0]; 2234d615ea0SBarry Smith last = 0; 2244d615ea0SBarry Smith for (i=1; i<=right; i++) { 2254d615ea0SBarry Smith if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);} 2264d615ea0SBarry Smith } 2274d615ea0SBarry Smith SWAP2(v[0],v[last],V[0],V[last],tmp); 2284d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v,V,last-1);CHKERRQ(ierr); 2294d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 2304d615ea0SBarry Smith PetscFunctionReturn(0); 2314d615ea0SBarry Smith } 2324d615ea0SBarry Smith 2334d615ea0SBarry Smith #undef __FUNCT__ 2344d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray" 2354d615ea0SBarry Smith /*@ 2364d615ea0SBarry Smith PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order; 2374d615ea0SBarry Smith changes a second array to match the sorted first array. 2384d615ea0SBarry Smith 2394d615ea0SBarry Smith Not Collective 2404d615ea0SBarry Smith 2414d615ea0SBarry Smith Input Parameters: 2424d615ea0SBarry Smith + n - number of values 2434d615ea0SBarry Smith . i - array of integers 2444d615ea0SBarry Smith - I - second array of integers 2454d615ea0SBarry Smith 2464d615ea0SBarry Smith Level: intermediate 2474d615ea0SBarry Smith 2484d615ea0SBarry Smith Concepts: sorting^ints with array 2494d615ea0SBarry Smith 2504d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt() 2514d615ea0SBarry Smith @*/ 2528738c821SJed Brown PetscErrorCode PETSCSYS_DLLEXPORT PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt i[],PetscMPIInt Ii[]) 2534d615ea0SBarry Smith { 2544d615ea0SBarry Smith PetscErrorCode ierr; 2554d615ea0SBarry Smith PetscMPIInt j,k,tmp,ik; 2564d615ea0SBarry Smith 2574d615ea0SBarry Smith PetscFunctionBegin; 2584d615ea0SBarry Smith if (n<8) { 2594d615ea0SBarry Smith for (k=0; k<n; k++) { 2604d615ea0SBarry Smith ik = i[k]; 2614d615ea0SBarry Smith for (j=k+1; j<n; j++) { 2624d615ea0SBarry Smith if (ik > i[j]) { 2634d615ea0SBarry Smith SWAP2(i[k],i[j],Ii[k],Ii[j],tmp); 2644d615ea0SBarry Smith ik = i[k]; 2654d615ea0SBarry Smith } 2664d615ea0SBarry Smith } 2674d615ea0SBarry Smith } 2684d615ea0SBarry Smith } else { 2694d615ea0SBarry Smith ierr = PetscSortMPIIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr); 2704d615ea0SBarry Smith } 2714d615ea0SBarry Smith PetscFunctionReturn(0); 2724d615ea0SBarry Smith } 2734d615ea0SBarry Smith 274e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/ 275e5c89e4eSSatish Balay #define SWAP2IntScalar(a,b,c,d,t,ts) {t=a;a=b;b=t;ts=c;c=d;d=ts;} 276e5c89e4eSSatish Balay 277e5c89e4eSSatish Balay #undef __FUNCT__ 278e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray_Private" 279e5c89e4eSSatish Balay /* 280e5c89e4eSSatish Balay Modified from PetscSortIntWithArray_Private(). 281e5c89e4eSSatish Balay */ 282e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithScalarArray_Private(PetscInt *v,PetscScalar *V,PetscInt right) 283e5c89e4eSSatish Balay { 284e5c89e4eSSatish Balay PetscErrorCode ierr; 285e5c89e4eSSatish Balay PetscInt i,vl,last,tmp; 286e5c89e4eSSatish Balay PetscScalar stmp; 287e5c89e4eSSatish Balay 288e5c89e4eSSatish Balay PetscFunctionBegin; 289e5c89e4eSSatish Balay if (right <= 1) { 290e5c89e4eSSatish Balay if (right == 1) { 291e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP2IntScalar(v[0],v[1],V[0],V[1],tmp,stmp); 292e5c89e4eSSatish Balay } 293e5c89e4eSSatish Balay PetscFunctionReturn(0); 294e5c89e4eSSatish Balay } 295e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[right/2],V[0],V[right/2],tmp,stmp); 296e5c89e4eSSatish Balay vl = v[0]; 297e5c89e4eSSatish Balay last = 0; 298e5c89e4eSSatish Balay for (i=1; i<=right; i++) { 299e5c89e4eSSatish Balay if (v[i] < vl) {last++; SWAP2IntScalar(v[last],v[i],V[last],V[i],tmp,stmp);} 300e5c89e4eSSatish Balay } 301e5c89e4eSSatish Balay SWAP2IntScalar(v[0],v[last],V[0],V[last],tmp,stmp); 302e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v,V,last-1);CHKERRQ(ierr); 303e5c89e4eSSatish Balay ierr = PetscSortIntWithScalarArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr); 304e5c89e4eSSatish Balay PetscFunctionReturn(0); 305e5c89e4eSSatish Balay } 306e5c89e4eSSatish Balay 307e5c89e4eSSatish Balay #undef __FUNCT__ 308e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray" 309e5c89e4eSSatish Balay /*@ 310e5c89e4eSSatish Balay PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order; 311e5c89e4eSSatish Balay changes a second SCALAR array to match the sorted first INTEGER array. 312e5c89e4eSSatish Balay 313e5c89e4eSSatish Balay Not Collective 314e5c89e4eSSatish Balay 315e5c89e4eSSatish Balay Input Parameters: 316e5c89e4eSSatish Balay + n - number of values 317e5c89e4eSSatish Balay . i - array of integers 318e5c89e4eSSatish Balay - I - second array of scalars 319e5c89e4eSSatish Balay 320e5c89e4eSSatish Balay Level: intermediate 321e5c89e4eSSatish Balay 322e5c89e4eSSatish Balay Concepts: sorting^ints with array 323e5c89e4eSSatish Balay 324e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray() 325e5c89e4eSSatish Balay @*/ 3268738c821SJed Brown PetscErrorCode PETSCSYS_DLLEXPORT PetscSortIntWithScalarArray(PetscInt n,PetscInt i[],PetscScalar Ii[]) 327e5c89e4eSSatish Balay { 328e5c89e4eSSatish Balay PetscErrorCode ierr; 329e5c89e4eSSatish Balay PetscInt j,k,tmp,ik; 330e5c89e4eSSatish Balay PetscScalar stmp; 331e5c89e4eSSatish Balay 332e5c89e4eSSatish Balay PetscFunctionBegin; 333e5c89e4eSSatish Balay if (n<8) { 334e5c89e4eSSatish Balay for (k=0; k<n; k++) { 335e5c89e4eSSatish Balay ik = i[k]; 336e5c89e4eSSatish Balay for (j=k+1; j<n; j++) { 337e5c89e4eSSatish Balay if (ik > i[j]) { 338b7940d39SSatish Balay SWAP2IntScalar(i[k],i[j],Ii[k],Ii[j],tmp,stmp); 339e5c89e4eSSatish Balay ik = i[k]; 340e5c89e4eSSatish Balay } 341e5c89e4eSSatish Balay } 342e5c89e4eSSatish Balay } 343e5c89e4eSSatish Balay } else { 344b7940d39SSatish Balay ierr = PetscSortIntWithScalarArray_Private(i,Ii,n-1);CHKERRQ(ierr); 345e5c89e4eSSatish Balay } 346e5c89e4eSSatish Balay PetscFunctionReturn(0); 347e5c89e4eSSatish Balay } 348e5c89e4eSSatish Balay 349e5c89e4eSSatish Balay 350*42eaa7f3SBarry Smith #undef __FUNCT__ 351*42eaa7f3SBarry Smith #define __FUNCT__ "PetscProcessTree" 352*42eaa7f3SBarry Smith /*@ 353*42eaa7f3SBarry Smith PetscProcessTree - Prepares tree data to be displayed graphically 354*42eaa7f3SBarry Smith 355*42eaa7f3SBarry Smith Not Collective 356*42eaa7f3SBarry Smith 357*42eaa7f3SBarry Smith Input Parameters: 358*42eaa7f3SBarry Smith + n - number of values 359*42eaa7f3SBarry Smith . mask - indicates those entries in the tree, location 0 is always masked 360*42eaa7f3SBarry Smith - parentid - indicates the parent of each entry 361*42eaa7f3SBarry Smith 362*42eaa7f3SBarry Smith Output Parameters: 363*42eaa7f3SBarry Smith + Nlevels - the number of levels 364*42eaa7f3SBarry Smith . Level - for each node tells its level 365*42eaa7f3SBarry Smith . Levelcnts - the number of nodes on each level 366*42eaa7f3SBarry Smith . Idbylevel - a list of ids on each of the levels, first level followed by second etc 367*42eaa7f3SBarry Smith - Column - for each id tells its column index 368*42eaa7f3SBarry Smith 369*42eaa7f3SBarry Smith Level: intermediate 370*42eaa7f3SBarry Smith 371*42eaa7f3SBarry Smith 372*42eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation() 373*42eaa7f3SBarry Smith @*/ 374*42eaa7f3SBarry Smith PetscErrorCode PETSCSYS_DLLEXPORT PetscProcessTree(PetscInt n,const PetscTruth mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column) 375*42eaa7f3SBarry Smith { 376*42eaa7f3SBarry Smith PetscInt i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column; 377*42eaa7f3SBarry Smith PetscErrorCode ierr; 378*42eaa7f3SBarry Smith PetscTruth done = PETSC_FALSE; 379*42eaa7f3SBarry Smith 380*42eaa7f3SBarry Smith PetscFunctionBegin; 381*42eaa7f3SBarry Smith if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set"); 382*42eaa7f3SBarry Smith for (i=0; i<n; i++) { 383*42eaa7f3SBarry Smith if (mask[i]) continue; 384*42eaa7f3SBarry Smith if (parentid[i] == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent"); 385*42eaa7f3SBarry Smith if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked"); 386*42eaa7f3SBarry Smith } 387*42eaa7f3SBarry Smith 388*42eaa7f3SBarry Smith 389*42eaa7f3SBarry Smith for (i=0; i<n; i++) { 390*42eaa7f3SBarry Smith if (!mask[i]) nmask++; 391*42eaa7f3SBarry Smith } 392*42eaa7f3SBarry Smith 393*42eaa7f3SBarry Smith /* determine the level in the tree of each node */ 394*42eaa7f3SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&level);CHKERRQ(ierr); 395*42eaa7f3SBarry Smith ierr = PetscMemzero(level,n*sizeof(PetscInt));CHKERRQ(ierr); 396*42eaa7f3SBarry Smith level[0] = 1; 397*42eaa7f3SBarry Smith while (!done) { 398*42eaa7f3SBarry Smith done = PETSC_TRUE; 399*42eaa7f3SBarry Smith for (i=0; i<n; i++) { 400*42eaa7f3SBarry Smith if (mask[i]) continue; 401*42eaa7f3SBarry Smith if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1; 402*42eaa7f3SBarry Smith else if (!level[i]) done = PETSC_FALSE; 403*42eaa7f3SBarry Smith } 404*42eaa7f3SBarry Smith } 405*42eaa7f3SBarry Smith for (i=0; i<n; i++) { 406*42eaa7f3SBarry Smith level[i]--; 407*42eaa7f3SBarry Smith nlevels = PetscMax(nlevels,level[i]); 408*42eaa7f3SBarry Smith } 409*42eaa7f3SBarry Smith 410*42eaa7f3SBarry Smith /* count the number of nodes on each level and its max */ 411*42eaa7f3SBarry Smith ierr = PetscMalloc(nlevels*sizeof(PetscInt),&levelcnt);CHKERRQ(ierr); 412*42eaa7f3SBarry Smith ierr = PetscMemzero(levelcnt,nlevels*sizeof(PetscInt));CHKERRQ(ierr); 413*42eaa7f3SBarry Smith for (i=0; i<n; i++) { 414*42eaa7f3SBarry Smith if (mask[i]) continue; 415*42eaa7f3SBarry Smith levelcnt[level[i]-1]++; 416*42eaa7f3SBarry Smith } 417*42eaa7f3SBarry Smith for (i=0; i<nlevels;i++) { 418*42eaa7f3SBarry Smith levelmax = PetscMax(levelmax,levelcnt[i]); 419*42eaa7f3SBarry Smith } 420*42eaa7f3SBarry Smith 421*42eaa7f3SBarry Smith /* for each level sort the ids by the parent id */ 422*42eaa7f3SBarry Smith ierr = PetscMalloc2(levelmax,PetscInt,&workid,levelmax,PetscInt,&workparentid);CHKERRQ(ierr); 423*42eaa7f3SBarry Smith ierr = PetscMalloc(nmask*sizeof(PetscInt),&idbylevel);CHKERRQ(ierr); 424*42eaa7f3SBarry Smith for (j=1; j<=nlevels;j++) { 425*42eaa7f3SBarry Smith cnt = 0; 426*42eaa7f3SBarry Smith for (i=0; i<n; i++) { 427*42eaa7f3SBarry Smith if (mask[i]) continue; 428*42eaa7f3SBarry Smith if (level[i] != j) continue; 429*42eaa7f3SBarry Smith workid[cnt] = i; 430*42eaa7f3SBarry Smith workparentid[cnt++] = parentid[i]; 431*42eaa7f3SBarry Smith } 432*42eaa7f3SBarry Smith /* PetscIntView(cnt,workparentid,0); 433*42eaa7f3SBarry Smith PetscIntView(cnt,workid,0); 434*42eaa7f3SBarry Smith ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr); 435*42eaa7f3SBarry Smith PetscIntView(cnt,workparentid,0); 436*42eaa7f3SBarry Smith PetscIntView(cnt,workid,0);*/ 437*42eaa7f3SBarry Smith ierr = PetscMemcpy(idbylevel+tcnt,workid,cnt*sizeof(PetscInt));CHKERRQ(ierr); 438*42eaa7f3SBarry Smith tcnt += cnt; 439*42eaa7f3SBarry Smith } 440*42eaa7f3SBarry Smith if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes"); 441*42eaa7f3SBarry Smith ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr); 442*42eaa7f3SBarry Smith 443*42eaa7f3SBarry Smith /* for each node list its column */ 444*42eaa7f3SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&column);CHKERRQ(ierr); 445*42eaa7f3SBarry Smith cnt = 0; 446*42eaa7f3SBarry Smith for (j=0; j<nlevels; j++) { 447*42eaa7f3SBarry Smith for (i=0; i<levelcnt[j]; i++) { 448*42eaa7f3SBarry Smith column[idbylevel[cnt++]] = i; 449*42eaa7f3SBarry Smith } 450*42eaa7f3SBarry Smith } 451*42eaa7f3SBarry Smith 452*42eaa7f3SBarry Smith *Nlevels = nlevels; 453*42eaa7f3SBarry Smith *Level = level; 454*42eaa7f3SBarry Smith *Levelcnt = levelcnt; 455*42eaa7f3SBarry Smith *Idbylevel = idbylevel; 456*42eaa7f3SBarry Smith *Column = column; 457*42eaa7f3SBarry Smith PetscFunctionReturn(0); 458*42eaa7f3SBarry Smith } 459