xref: /petsc/src/sys/utils/sorti.c (revision d9f1162d576517911b479a3b323db8c5eb0c1109)
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 /*@
129*d9f1162dSJed Brown   PetscFindInt - Finds integer 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
135*d9f1162dSJed Brown .  n   - number of values in the array
13660e03357SMatthew G Knepley -  ii  - array of integers
13760e03357SMatthew G Knepley 
13860e03357SMatthew G Knepley    Output Parameter:
139*d9f1162dSJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
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 @*/
147026ec6cbSMatthew G Knepley PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt ii[], PetscInt *loc)
14860e03357SMatthew G Knepley {
149*d9f1162dSJed Brown   PetscInt lo = 0,hi = n;
15060e03357SMatthew G Knepley 
15160e03357SMatthew G Knepley   PetscFunctionBegin;
15260e03357SMatthew G Knepley   PetscValidPointer(ii,3);
15360e03357SMatthew G Knepley   PetscValidPointer(loc,4);
154*d9f1162dSJed Brown   while (hi - lo > 1) {
155*d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo)/2;
156*d9f1162dSJed Brown     if (key < ii[mid]) hi = mid;
157*d9f1162dSJed Brown     else               lo = mid;
15860e03357SMatthew G Knepley   }
159*d9f1162dSJed Brown   *loc = key == ii[lo] ? lo : -(lo + (key > ii[lo]) + 1);
16060e03357SMatthew G Knepley   PetscFunctionReturn(0);
16160e03357SMatthew G Knepley }
16260e03357SMatthew G Knepley 
1631db36a52SBarry Smith 
164e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/
165e5c89e4eSSatish Balay #define SWAP2(a,b,c,d,t) {t=a;a=b;b=t;t=c;c=d;d=t;}
166e5c89e4eSSatish Balay 
167e5c89e4eSSatish Balay #undef __FUNCT__
168e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithArray_Private"
169e5c89e4eSSatish Balay /*
170e5c89e4eSSatish Balay    A simple version of quicksort; taken from Kernighan and Ritchie, page 87.
171e5c89e4eSSatish Balay    Assumes 0 origin for v, number of elements = right+1 (right is index of
172e5c89e4eSSatish Balay    right-most member).
173e5c89e4eSSatish Balay */
174e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithArray_Private(PetscInt *v,PetscInt *V,PetscInt right)
175e5c89e4eSSatish Balay {
176e5c89e4eSSatish Balay   PetscErrorCode ierr;
177e5c89e4eSSatish Balay   PetscInt       i,vl,last,tmp;
178e5c89e4eSSatish Balay 
179e5c89e4eSSatish Balay   PetscFunctionBegin;
180e5c89e4eSSatish Balay   if (right <= 1) {
181e5c89e4eSSatish Balay     if (right == 1) {
182e5c89e4eSSatish Balay       if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp);
183e5c89e4eSSatish Balay     }
184e5c89e4eSSatish Balay     PetscFunctionReturn(0);
185e5c89e4eSSatish Balay   }
186e5c89e4eSSatish Balay   SWAP2(v[0],v[right/2],V[0],V[right/2],tmp);
187e5c89e4eSSatish Balay   vl   = v[0];
188e5c89e4eSSatish Balay   last = 0;
189e5c89e4eSSatish Balay   for (i=1; i<=right; i++) {
190e5c89e4eSSatish Balay     if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);}
191e5c89e4eSSatish Balay   }
192e5c89e4eSSatish Balay   SWAP2(v[0],v[last],V[0],V[last],tmp);
193e5c89e4eSSatish Balay   ierr = PetscSortIntWithArray_Private(v,V,last-1);CHKERRQ(ierr);
194e5c89e4eSSatish Balay   ierr = PetscSortIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr);
195e5c89e4eSSatish Balay   PetscFunctionReturn(0);
196e5c89e4eSSatish Balay }
197e5c89e4eSSatish Balay 
198e5c89e4eSSatish Balay #undef __FUNCT__
199e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithArray"
200e5c89e4eSSatish Balay /*@
201e5c89e4eSSatish Balay    PetscSortIntWithArray - Sorts an array of integers in place in increasing order;
202e5c89e4eSSatish Balay        changes a second array to match the sorted first array.
203e5c89e4eSSatish Balay 
204e5c89e4eSSatish Balay    Not Collective
205e5c89e4eSSatish Balay 
206e5c89e4eSSatish Balay    Input Parameters:
207e5c89e4eSSatish Balay +  n  - number of values
208e5c89e4eSSatish Balay .  i  - array of integers
209e5c89e4eSSatish Balay -  I - second array of integers
210e5c89e4eSSatish Balay 
211e5c89e4eSSatish Balay    Level: intermediate
212e5c89e4eSSatish Balay 
213e5c89e4eSSatish Balay    Concepts: sorting^ints with array
214e5c89e4eSSatish Balay 
215e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
216e5c89e4eSSatish Balay @*/
2177087cfbeSBarry Smith PetscErrorCode  PetscSortIntWithArray(PetscInt n,PetscInt i[],PetscInt Ii[])
218e5c89e4eSSatish Balay {
219e5c89e4eSSatish Balay   PetscErrorCode ierr;
220e5c89e4eSSatish Balay   PetscInt       j,k,tmp,ik;
221e5c89e4eSSatish Balay 
222e5c89e4eSSatish Balay   PetscFunctionBegin;
223e5c89e4eSSatish Balay   if (n<8) {
224e5c89e4eSSatish Balay     for (k=0; k<n; k++) {
225e5c89e4eSSatish Balay       ik = i[k];
226e5c89e4eSSatish Balay       for (j=k+1; j<n; j++) {
227e5c89e4eSSatish Balay 	if (ik > i[j]) {
228b7940d39SSatish Balay 	  SWAP2(i[k],i[j],Ii[k],Ii[j],tmp);
229e5c89e4eSSatish Balay 	  ik = i[k];
230e5c89e4eSSatish Balay 	}
231e5c89e4eSSatish Balay       }
232e5c89e4eSSatish Balay     }
233e5c89e4eSSatish Balay   } else {
234b7940d39SSatish Balay     ierr = PetscSortIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr);
235e5c89e4eSSatish Balay   }
236e5c89e4eSSatish Balay   PetscFunctionReturn(0);
237e5c89e4eSSatish Balay }
238e5c89e4eSSatish Balay 
239c1f0200aSDmitry Karpeev 
240c1f0200aSDmitry 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;}
241c1f0200aSDmitry Karpeev 
242c1f0200aSDmitry Karpeev #undef __FUNCT__
243c1f0200aSDmitry Karpeev #define __FUNCT__ "PetscSortIntWithArrayPair_Private"
244c1f0200aSDmitry Karpeev /*
245c1f0200aSDmitry Karpeev    A simple version of quicksort; taken from Kernighan and Ritchie, page 87.
246c1f0200aSDmitry Karpeev    Assumes 0 origin for v, number of elements = right+1 (right is index of
247c1f0200aSDmitry Karpeev    right-most member).
248c1f0200aSDmitry Karpeev */
249d7aa01a8SHong Zhang static PetscErrorCode PetscSortIntWithArrayPair_Private(PetscInt *L,PetscInt *J, PetscInt *K, PetscInt right)
250c1f0200aSDmitry Karpeev {
251c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
252c1f0200aSDmitry Karpeev   PetscInt       i,vl,last,tmp;
253c1f0200aSDmitry Karpeev 
254c1f0200aSDmitry Karpeev   PetscFunctionBegin;
255c1f0200aSDmitry Karpeev   if (right <= 1) {
256c1f0200aSDmitry Karpeev     if (right == 1) {
257d7aa01a8SHong Zhang       if (L[0] > L[1]) SWAP3(L[0],L[1],J[0],J[1],K[0],K[1],tmp);
258c1f0200aSDmitry Karpeev     }
259c1f0200aSDmitry Karpeev     PetscFunctionReturn(0);
260c1f0200aSDmitry Karpeev   }
261d7aa01a8SHong Zhang   SWAP3(L[0],L[right/2],J[0],J[right/2],K[0],K[right/2],tmp);
262d7aa01a8SHong Zhang   vl   = L[0];
263c1f0200aSDmitry Karpeev   last = 0;
264c1f0200aSDmitry Karpeev   for (i=1; i<=right; i++) {
265d7aa01a8SHong Zhang     if (L[i] < vl) {last++; SWAP3(L[last],L[i],J[last],J[i],K[last],K[i],tmp);}
266c1f0200aSDmitry Karpeev   }
267d7aa01a8SHong Zhang   SWAP3(L[0],L[last],J[0],J[last],K[0],K[last],tmp);
268d7aa01a8SHong Zhang   ierr = PetscSortIntWithArrayPair_Private(L,J,K,last-1);CHKERRQ(ierr);
269d7aa01a8SHong Zhang   ierr = PetscSortIntWithArrayPair_Private(L+last+1,J+last+1,K+last+1,right-(last+1));CHKERRQ(ierr);
270c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
271c1f0200aSDmitry Karpeev }
272c1f0200aSDmitry Karpeev 
273c1f0200aSDmitry Karpeev #undef __FUNCT__
274c1f0200aSDmitry Karpeev #define __FUNCT__ "PetscSortIntWithArrayPair"
275c1f0200aSDmitry Karpeev /*@
276c1f0200aSDmitry Karpeev    PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order;
277c1f0200aSDmitry Karpeev        changes a pair of integer arrays to match the sorted first array.
278c1f0200aSDmitry Karpeev 
279c1f0200aSDmitry Karpeev    Not Collective
280c1f0200aSDmitry Karpeev 
281c1f0200aSDmitry Karpeev    Input Parameters:
282c1f0200aSDmitry Karpeev +  n  - number of values
283c1f0200aSDmitry Karpeev .  I  - array of integers
284c1f0200aSDmitry Karpeev .  J  - second array of integers (first array of the pair)
285c1f0200aSDmitry Karpeev -  K  - third array of integers  (second array of the pair)
286c1f0200aSDmitry Karpeev 
287c1f0200aSDmitry Karpeev    Level: intermediate
288c1f0200aSDmitry Karpeev 
289c1f0200aSDmitry Karpeev    Concepts: sorting^ints with array pair
290c1f0200aSDmitry Karpeev 
291c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray()
292c1f0200aSDmitry Karpeev @*/
293d7aa01a8SHong Zhang PetscErrorCode  PetscSortIntWithArrayPair(PetscInt n,PetscInt *L,PetscInt *J, PetscInt *K)
294c1f0200aSDmitry Karpeev {
295c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
296c1f0200aSDmitry Karpeev   PetscInt       j,k,tmp,ik;
297c1f0200aSDmitry Karpeev 
298c1f0200aSDmitry Karpeev   PetscFunctionBegin;
299c1f0200aSDmitry Karpeev   if (n<8) {
300c1f0200aSDmitry Karpeev     for (k=0; k<n; k++) {
301d7aa01a8SHong Zhang       ik = L[k];
302c1f0200aSDmitry Karpeev       for (j=k+1; j<n; j++) {
303d7aa01a8SHong Zhang 	if (ik > L[j]) {
304d7aa01a8SHong Zhang 	  SWAP3(L[k],L[j],J[k],J[j],K[k],K[j],tmp);
305d7aa01a8SHong Zhang 	  ik = L[k];
306c1f0200aSDmitry Karpeev 	}
307c1f0200aSDmitry Karpeev       }
308c1f0200aSDmitry Karpeev     }
309c1f0200aSDmitry Karpeev   } else {
310d7aa01a8SHong Zhang     ierr = PetscSortIntWithArrayPair_Private(L,J,K,n-1);CHKERRQ(ierr);
311c1f0200aSDmitry Karpeev   }
312c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
313c1f0200aSDmitry Karpeev }
314c1f0200aSDmitry Karpeev 
31517d7d925SStefano Zampini #undef __FUNCT__
31617d7d925SStefano Zampini #define __FUNCT__ "PetscSortMPIInt_Private"
31717d7d925SStefano Zampini /*
31817d7d925SStefano Zampini    A simple version of quicksort; taken from Kernighan and Ritchie, page 87.
31917d7d925SStefano Zampini    Assumes 0 origin for v, number of elements = right+1 (right is index of
32017d7d925SStefano Zampini    right-most member).
32117d7d925SStefano Zampini */
32217d7d925SStefano Zampini static void PetscSortMPIInt_Private(PetscMPIInt *v,PetscInt right)
32317d7d925SStefano Zampini {
32417d7d925SStefano Zampini   PetscInt          i,j;
32517d7d925SStefano Zampini   PetscMPIInt       pivot,tmp;
32617d7d925SStefano Zampini 
32717d7d925SStefano Zampini   if (right <= 1) {
32817d7d925SStefano Zampini     if (right == 1) {
32917d7d925SStefano Zampini       if (v[0] > v[1]) SWAP(v[0],v[1],tmp);
33017d7d925SStefano Zampini     }
33117d7d925SStefano Zampini     return;
33217d7d925SStefano Zampini   }
33317d7d925SStefano Zampini   i = MEDIAN(v,right);          /* Choose a pivot */
33417d7d925SStefano Zampini   SWAP(v[0],v[i],tmp);          /* Move it out of the way */
33517d7d925SStefano Zampini   pivot = v[0];
33617d7d925SStefano Zampini   for (i=0,j=right+1;;) {
33717d7d925SStefano Zampini     while (++i < j && v[i] <= pivot); /* Scan from the left */
33817d7d925SStefano Zampini     while (v[--j] > pivot);           /* Scan from the right */
33917d7d925SStefano Zampini     if (i >= j) break;
34017d7d925SStefano Zampini     SWAP(v[i],v[j],tmp);
34117d7d925SStefano Zampini   }
34217d7d925SStefano Zampini   SWAP(v[0],v[j],tmp);          /* Put pivot back in place. */
34317d7d925SStefano Zampini   PetscSortMPIInt_Private(v,j-1);
34417d7d925SStefano Zampini   PetscSortMPIInt_Private(v+j+1,right-(j+1));
34517d7d925SStefano Zampini }
34617d7d925SStefano Zampini 
34717d7d925SStefano Zampini #undef __FUNCT__
34817d7d925SStefano Zampini #define __FUNCT__ "PetscSortMPIInt"
34917d7d925SStefano Zampini /*@
35017d7d925SStefano Zampini    PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order.
35117d7d925SStefano Zampini 
35217d7d925SStefano Zampini    Not Collective
35317d7d925SStefano Zampini 
35417d7d925SStefano Zampini    Input Parameters:
35517d7d925SStefano Zampini +  n  - number of values
35617d7d925SStefano Zampini -  i  - array of integers
35717d7d925SStefano Zampini 
35817d7d925SStefano Zampini    Level: intermediate
35917d7d925SStefano Zampini 
36017d7d925SStefano Zampini    Concepts: sorting^ints
36117d7d925SStefano Zampini 
36217d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation()
36317d7d925SStefano Zampini @*/
36417d7d925SStefano Zampini PetscErrorCode  PetscSortMPIInt(PetscInt n,PetscMPIInt i[])
36517d7d925SStefano Zampini {
36617d7d925SStefano Zampini   PetscInt       j,k;
36717d7d925SStefano Zampini   PetscMPIInt    tmp,ik;
36817d7d925SStefano Zampini 
36917d7d925SStefano Zampini   PetscFunctionBegin;
37017d7d925SStefano Zampini   if (n<8) {
37117d7d925SStefano Zampini     for (k=0; k<n; k++) {
37217d7d925SStefano Zampini       ik = i[k];
37317d7d925SStefano Zampini       for (j=k+1; j<n; j++) {
37417d7d925SStefano Zampini 	if (ik > i[j]) {
37517d7d925SStefano Zampini 	  SWAP(i[k],i[j],tmp);
37617d7d925SStefano Zampini 	  ik = i[k];
37717d7d925SStefano Zampini 	}
37817d7d925SStefano Zampini       }
37917d7d925SStefano Zampini     }
38017d7d925SStefano Zampini   } else {
38117d7d925SStefano Zampini     PetscSortMPIInt_Private(i,n-1);
38217d7d925SStefano Zampini   }
38317d7d925SStefano Zampini   PetscFunctionReturn(0);
38417d7d925SStefano Zampini }
38517d7d925SStefano Zampini 
38617d7d925SStefano Zampini #undef __FUNCT__
38717d7d925SStefano Zampini #define __FUNCT__ "PetscSortRemoveDupsMPIInt"
38817d7d925SStefano Zampini /*@
38917d7d925SStefano Zampini    PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries
39017d7d925SStefano Zampini 
39117d7d925SStefano Zampini    Not Collective
39217d7d925SStefano Zampini 
39317d7d925SStefano Zampini    Input Parameters:
39417d7d925SStefano Zampini +  n  - number of values
39517d7d925SStefano Zampini -  ii  - array of integers
39617d7d925SStefano Zampini 
39717d7d925SStefano Zampini    Output Parameter:
39817d7d925SStefano Zampini .  n - number of non-redundant values
39917d7d925SStefano Zampini 
40017d7d925SStefano Zampini    Level: intermediate
40117d7d925SStefano Zampini 
40217d7d925SStefano Zampini    Concepts: sorting^ints
40317d7d925SStefano Zampini 
40417d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
40517d7d925SStefano Zampini @*/
40617d7d925SStefano Zampini PetscErrorCode  PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt ii[])
40717d7d925SStefano Zampini {
40817d7d925SStefano Zampini   PetscErrorCode ierr;
40917d7d925SStefano Zampini   PetscInt       i,s = 0,N = *n, b = 0;
41017d7d925SStefano Zampini 
41117d7d925SStefano Zampini   PetscFunctionBegin;
41217d7d925SStefano Zampini   ierr = PetscSortMPIInt(N,ii);CHKERRQ(ierr);
41317d7d925SStefano Zampini   for (i=0; i<N-1; i++) {
41417d7d925SStefano Zampini     if (ii[b+s+1] != ii[b]) {ii[b+1] = ii[b+s+1]; b++;}
41517d7d925SStefano Zampini     else s++;
41617d7d925SStefano Zampini   }
41717d7d925SStefano Zampini   *n = N - s;
41817d7d925SStefano Zampini   PetscFunctionReturn(0);
41917d7d925SStefano Zampini }
420c1f0200aSDmitry Karpeev 
4214d615ea0SBarry Smith #undef __FUNCT__
4224d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray_Private"
4234d615ea0SBarry Smith /*
4244d615ea0SBarry Smith    A simple version of quicksort; taken from Kernighan and Ritchie, page 87.
4254d615ea0SBarry Smith    Assumes 0 origin for v, number of elements = right+1 (right is index of
4264d615ea0SBarry Smith    right-most member).
4274d615ea0SBarry Smith */
4284d615ea0SBarry Smith static PetscErrorCode PetscSortMPIIntWithArray_Private(PetscMPIInt *v,PetscMPIInt *V,PetscMPIInt right)
4294d615ea0SBarry Smith {
4304d615ea0SBarry Smith   PetscErrorCode ierr;
4314d615ea0SBarry Smith   PetscMPIInt    i,vl,last,tmp;
4324d615ea0SBarry Smith 
4334d615ea0SBarry Smith   PetscFunctionBegin;
4344d615ea0SBarry Smith   if (right <= 1) {
4354d615ea0SBarry Smith     if (right == 1) {
4364d615ea0SBarry Smith       if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp);
4374d615ea0SBarry Smith     }
4384d615ea0SBarry Smith     PetscFunctionReturn(0);
4394d615ea0SBarry Smith   }
4404d615ea0SBarry Smith   SWAP2(v[0],v[right/2],V[0],V[right/2],tmp);
4414d615ea0SBarry Smith   vl   = v[0];
4424d615ea0SBarry Smith   last = 0;
4434d615ea0SBarry Smith   for (i=1; i<=right; i++) {
4444d615ea0SBarry Smith     if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);}
4454d615ea0SBarry Smith   }
4464d615ea0SBarry Smith   SWAP2(v[0],v[last],V[0],V[last],tmp);
4474d615ea0SBarry Smith   ierr = PetscSortMPIIntWithArray_Private(v,V,last-1);CHKERRQ(ierr);
4484d615ea0SBarry Smith   ierr = PetscSortMPIIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr);
4494d615ea0SBarry Smith   PetscFunctionReturn(0);
4504d615ea0SBarry Smith }
4514d615ea0SBarry Smith 
4524d615ea0SBarry Smith #undef __FUNCT__
4534d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray"
4544d615ea0SBarry Smith /*@
4554d615ea0SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order;
4564d615ea0SBarry Smith        changes a second array to match the sorted first array.
4574d615ea0SBarry Smith 
4584d615ea0SBarry Smith    Not Collective
4594d615ea0SBarry Smith 
4604d615ea0SBarry Smith    Input Parameters:
4614d615ea0SBarry Smith +  n  - number of values
4624d615ea0SBarry Smith .  i  - array of integers
4634d615ea0SBarry Smith -  I - second array of integers
4644d615ea0SBarry Smith 
4654d615ea0SBarry Smith    Level: intermediate
4664d615ea0SBarry Smith 
4674d615ea0SBarry Smith    Concepts: sorting^ints with array
4684d615ea0SBarry Smith 
4694d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
4704d615ea0SBarry Smith @*/
4717087cfbeSBarry Smith PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt i[],PetscMPIInt Ii[])
4724d615ea0SBarry Smith {
4734d615ea0SBarry Smith   PetscErrorCode ierr;
4744d615ea0SBarry Smith   PetscMPIInt    j,k,tmp,ik;
4754d615ea0SBarry Smith 
4764d615ea0SBarry Smith   PetscFunctionBegin;
4774d615ea0SBarry Smith   if (n<8) {
4784d615ea0SBarry Smith     for (k=0; k<n; k++) {
4794d615ea0SBarry Smith       ik = i[k];
4804d615ea0SBarry Smith       for (j=k+1; j<n; j++) {
4814d615ea0SBarry Smith 	if (ik > i[j]) {
4824d615ea0SBarry Smith 	  SWAP2(i[k],i[j],Ii[k],Ii[j],tmp);
4834d615ea0SBarry Smith 	  ik = i[k];
4844d615ea0SBarry Smith 	}
4854d615ea0SBarry Smith       }
4864d615ea0SBarry Smith     }
4874d615ea0SBarry Smith   } else {
4884d615ea0SBarry Smith     ierr = PetscSortMPIIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr);
4894d615ea0SBarry Smith   }
4904d615ea0SBarry Smith   PetscFunctionReturn(0);
4914d615ea0SBarry Smith }
4924d615ea0SBarry Smith 
493e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/
494e5c89e4eSSatish Balay #define SWAP2IntScalar(a,b,c,d,t,ts) {t=a;a=b;b=t;ts=c;c=d;d=ts;}
495e5c89e4eSSatish Balay 
496e5c89e4eSSatish Balay #undef __FUNCT__
497e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray_Private"
498e5c89e4eSSatish Balay /*
499e5c89e4eSSatish Balay    Modified from PetscSortIntWithArray_Private().
500e5c89e4eSSatish Balay */
501e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithScalarArray_Private(PetscInt *v,PetscScalar *V,PetscInt right)
502e5c89e4eSSatish Balay {
503e5c89e4eSSatish Balay   PetscErrorCode ierr;
504e5c89e4eSSatish Balay   PetscInt       i,vl,last,tmp;
505e5c89e4eSSatish Balay   PetscScalar    stmp;
506e5c89e4eSSatish Balay 
507e5c89e4eSSatish Balay   PetscFunctionBegin;
508e5c89e4eSSatish Balay   if (right <= 1) {
509e5c89e4eSSatish Balay     if (right == 1) {
510e5c89e4eSSatish Balay       if (v[0] > v[1]) SWAP2IntScalar(v[0],v[1],V[0],V[1],tmp,stmp);
511e5c89e4eSSatish Balay     }
512e5c89e4eSSatish Balay     PetscFunctionReturn(0);
513e5c89e4eSSatish Balay   }
514e5c89e4eSSatish Balay   SWAP2IntScalar(v[0],v[right/2],V[0],V[right/2],tmp,stmp);
515e5c89e4eSSatish Balay   vl   = v[0];
516e5c89e4eSSatish Balay   last = 0;
517e5c89e4eSSatish Balay   for (i=1; i<=right; i++) {
518e5c89e4eSSatish Balay     if (v[i] < vl) {last++; SWAP2IntScalar(v[last],v[i],V[last],V[i],tmp,stmp);}
519e5c89e4eSSatish Balay   }
520e5c89e4eSSatish Balay   SWAP2IntScalar(v[0],v[last],V[0],V[last],tmp,stmp);
521e5c89e4eSSatish Balay   ierr = PetscSortIntWithScalarArray_Private(v,V,last-1);CHKERRQ(ierr);
522e5c89e4eSSatish Balay   ierr = PetscSortIntWithScalarArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr);
523e5c89e4eSSatish Balay   PetscFunctionReturn(0);
524e5c89e4eSSatish Balay }
525e5c89e4eSSatish Balay 
526e5c89e4eSSatish Balay #undef __FUNCT__
527e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray"
528e5c89e4eSSatish Balay /*@
529e5c89e4eSSatish Balay    PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
530e5c89e4eSSatish Balay        changes a second SCALAR array to match the sorted first INTEGER array.
531e5c89e4eSSatish Balay 
532e5c89e4eSSatish Balay    Not Collective
533e5c89e4eSSatish Balay 
534e5c89e4eSSatish Balay    Input Parameters:
535e5c89e4eSSatish Balay +  n  - number of values
536e5c89e4eSSatish Balay .  i  - array of integers
537e5c89e4eSSatish Balay -  I - second array of scalars
538e5c89e4eSSatish Balay 
539e5c89e4eSSatish Balay    Level: intermediate
540e5c89e4eSSatish Balay 
541e5c89e4eSSatish Balay    Concepts: sorting^ints with array
542e5c89e4eSSatish Balay 
543e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
544e5c89e4eSSatish Balay @*/
5457087cfbeSBarry Smith PetscErrorCode  PetscSortIntWithScalarArray(PetscInt n,PetscInt i[],PetscScalar Ii[])
546e5c89e4eSSatish Balay {
547e5c89e4eSSatish Balay   PetscErrorCode ierr;
548e5c89e4eSSatish Balay   PetscInt       j,k,tmp,ik;
549e5c89e4eSSatish Balay   PetscScalar    stmp;
550e5c89e4eSSatish Balay 
551e5c89e4eSSatish Balay   PetscFunctionBegin;
552e5c89e4eSSatish Balay   if (n<8) {
553e5c89e4eSSatish Balay     for (k=0; k<n; k++) {
554e5c89e4eSSatish Balay       ik = i[k];
555e5c89e4eSSatish Balay       for (j=k+1; j<n; j++) {
556e5c89e4eSSatish Balay 	if (ik > i[j]) {
557b7940d39SSatish Balay 	  SWAP2IntScalar(i[k],i[j],Ii[k],Ii[j],tmp,stmp);
558e5c89e4eSSatish Balay 	  ik = i[k];
559e5c89e4eSSatish Balay 	}
560e5c89e4eSSatish Balay       }
561e5c89e4eSSatish Balay     }
562e5c89e4eSSatish Balay   } else {
563b7940d39SSatish Balay     ierr = PetscSortIntWithScalarArray_Private(i,Ii,n-1);CHKERRQ(ierr);
564e5c89e4eSSatish Balay   }
565e5c89e4eSSatish Balay   PetscFunctionReturn(0);
566e5c89e4eSSatish Balay }
567e5c89e4eSSatish Balay 
568b4301105SBarry Smith 
569b4301105SBarry Smith 
570b4301105SBarry Smith #undef __FUNCT__
571c1f0200aSDmitry Karpeev #define __FUNCT__ "PetscMergeIntArrayPair"
572c1f0200aSDmitry Karpeev /*@
573c1f0200aSDmitry Karpeev    PetscMergeIntArrayPair -     Merges two SORTED integer arrays along with an additional array of integers.
574c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
575c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
576c1f0200aSDmitry Karpeev 
577c1f0200aSDmitry Karpeev    Not Collective
578c1f0200aSDmitry Karpeev 
579c1f0200aSDmitry Karpeev    Input Parameters:
580c1f0200aSDmitry Karpeev +  an  - number of values in the first array
581c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
582c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
583c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
584c1f0200aSDmitry Karpeev .  bI  - second array of integers
585c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
586c1f0200aSDmitry Karpeev 
587c1f0200aSDmitry Karpeev    Output Parameters:
588c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
589c1f0200aSDmitry Karpeev .  I   - merged sorted array
590c1f0200aSDmitry Karpeev -  J   - merged additional array
591c1f0200aSDmitry Karpeev 
592c1f0200aSDmitry Karpeev    Level: intermediate
593c1f0200aSDmitry Karpeev 
594c1f0200aSDmitry Karpeev    Concepts: merging^arrays
595c1f0200aSDmitry Karpeev 
596c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
597c1f0200aSDmitry Karpeev @*/
598d7aa01a8SHong 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)
599c1f0200aSDmitry Karpeev {
600c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
601d7aa01a8SHong Zhang   PetscInt n_, *L_ = *L, *J_= *J, ak, bk, k;
602c1f0200aSDmitry Karpeev 
603c1f0200aSDmitry Karpeev   n_ = an + bn;
604c1f0200aSDmitry Karpeev   *n = n_;
605d7aa01a8SHong Zhang   if (!L_) {
606d7aa01a8SHong Zhang     ierr = PetscMalloc(n_*sizeof(PetscInt), L); CHKERRQ(ierr);
607d7aa01a8SHong Zhang     L_ = *L;
608c1f0200aSDmitry Karpeev   }
609c1f0200aSDmitry Karpeev   if (!J_){
610c1f0200aSDmitry Karpeev     ierr = PetscMalloc(n_*sizeof(PetscInt), &J_); CHKERRQ(ierr);
611c1f0200aSDmitry Karpeev     J_ = *J;
612c1f0200aSDmitry Karpeev   }
613c1f0200aSDmitry Karpeev   k = ak = bk = 0;
614c1f0200aSDmitry Karpeev   while(ak < an && bk < bn) {
615c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
616d7aa01a8SHong Zhang       L_[k] = aI[ak];
617c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
618c1f0200aSDmitry Karpeev       ++ak;
619c1f0200aSDmitry Karpeev       ++k;
620c1f0200aSDmitry Karpeev     }
621c1f0200aSDmitry Karpeev     else {
622d7aa01a8SHong Zhang       L_[k] = bI[bk];
623c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
624c1f0200aSDmitry Karpeev       ++bk;
625c1f0200aSDmitry Karpeev       ++k;
626c1f0200aSDmitry Karpeev     }
627c1f0200aSDmitry Karpeev   }
628c1f0200aSDmitry Karpeev   if (ak < an) {
629d7aa01a8SHong Zhang     ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt)); CHKERRQ(ierr);
630c1f0200aSDmitry Karpeev     ierr = PetscMemcpy(J_+k,aJ+ak,(an-ak)*sizeof(PetscInt)); CHKERRQ(ierr);
631c1f0200aSDmitry Karpeev     k += (an-ak);
632c1f0200aSDmitry Karpeev   }
633c1f0200aSDmitry Karpeev   if (bk < bn) {
634d7aa01a8SHong Zhang     ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt)); CHKERRQ(ierr);
635c1f0200aSDmitry Karpeev     ierr = PetscMemcpy(J_+k,bJ+bk,(bn-bk)*sizeof(PetscInt)); CHKERRQ(ierr);
636c1f0200aSDmitry Karpeev     k += (bn-bk);
637c1f0200aSDmitry Karpeev   }
638c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
639c1f0200aSDmitry Karpeev }
640c1f0200aSDmitry Karpeev 
641e5c89e4eSSatish Balay 
64242eaa7f3SBarry Smith #undef __FUNCT__
64342eaa7f3SBarry Smith #define __FUNCT__ "PetscProcessTree"
64442eaa7f3SBarry Smith /*@
64542eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
64642eaa7f3SBarry Smith 
64742eaa7f3SBarry Smith    Not Collective
64842eaa7f3SBarry Smith 
64942eaa7f3SBarry Smith    Input Parameters:
65042eaa7f3SBarry Smith +  n  - number of values
65142eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
65242eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
65342eaa7f3SBarry Smith 
65442eaa7f3SBarry Smith    Output Parameters:
65542eaa7f3SBarry Smith +  Nlevels - the number of levels
65642eaa7f3SBarry Smith .  Level - for each node tells its level
65742eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
65842eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
65942eaa7f3SBarry Smith -  Column - for each id tells its column index
66042eaa7f3SBarry Smith 
66142eaa7f3SBarry Smith    Level: intermediate
66242eaa7f3SBarry Smith 
66342eaa7f3SBarry Smith 
66442eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation()
66542eaa7f3SBarry Smith @*/
6667087cfbeSBarry Smith PetscErrorCode  PetscProcessTree(PetscInt n,const PetscBool  mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column)
66742eaa7f3SBarry Smith {
66842eaa7f3SBarry Smith   PetscInt       i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column;
66942eaa7f3SBarry Smith   PetscErrorCode ierr;
670ace3abfcSBarry Smith   PetscBool      done = PETSC_FALSE;
67142eaa7f3SBarry Smith 
67242eaa7f3SBarry Smith   PetscFunctionBegin;
67342eaa7f3SBarry Smith   if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set");
67442eaa7f3SBarry Smith   for (i=0; i<n; i++) {
67542eaa7f3SBarry Smith     if (mask[i]) continue;
67642eaa7f3SBarry Smith     if (parentid[i]  == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent");
67742eaa7f3SBarry Smith     if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked");
67842eaa7f3SBarry Smith   }
67942eaa7f3SBarry Smith 
68042eaa7f3SBarry Smith 
68142eaa7f3SBarry Smith   for (i=0; i<n; i++) {
68242eaa7f3SBarry Smith     if (!mask[i]) nmask++;
68342eaa7f3SBarry Smith   }
68442eaa7f3SBarry Smith 
68542eaa7f3SBarry Smith   /* determine the level in the tree of each node */
68642eaa7f3SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&level);CHKERRQ(ierr);
68742eaa7f3SBarry Smith   ierr = PetscMemzero(level,n*sizeof(PetscInt));CHKERRQ(ierr);
68842eaa7f3SBarry Smith   level[0] = 1;
68942eaa7f3SBarry Smith   while (!done) {
69042eaa7f3SBarry Smith     done = PETSC_TRUE;
69142eaa7f3SBarry Smith     for (i=0; i<n; i++) {
69242eaa7f3SBarry Smith       if (mask[i]) continue;
69342eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
69442eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
69542eaa7f3SBarry Smith     }
69642eaa7f3SBarry Smith   }
69742eaa7f3SBarry Smith   for (i=0; i<n; i++) {
69842eaa7f3SBarry Smith     level[i]--;
69942eaa7f3SBarry Smith     nlevels = PetscMax(nlevels,level[i]);
70042eaa7f3SBarry Smith   }
70142eaa7f3SBarry Smith 
70242eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
70342eaa7f3SBarry Smith   ierr = PetscMalloc(nlevels*sizeof(PetscInt),&levelcnt);CHKERRQ(ierr);
70442eaa7f3SBarry Smith   ierr = PetscMemzero(levelcnt,nlevels*sizeof(PetscInt));CHKERRQ(ierr);
70542eaa7f3SBarry Smith   for (i=0; i<n; i++) {
70642eaa7f3SBarry Smith     if (mask[i]) continue;
70742eaa7f3SBarry Smith     levelcnt[level[i]-1]++;
70842eaa7f3SBarry Smith   }
70942eaa7f3SBarry Smith   for (i=0; i<nlevels;i++) {
71042eaa7f3SBarry Smith     levelmax = PetscMax(levelmax,levelcnt[i]);
71142eaa7f3SBarry Smith   }
71242eaa7f3SBarry Smith 
71342eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
71442eaa7f3SBarry Smith   ierr = PetscMalloc2(levelmax,PetscInt,&workid,levelmax,PetscInt,&workparentid);CHKERRQ(ierr);
71542eaa7f3SBarry Smith   ierr = PetscMalloc(nmask*sizeof(PetscInt),&idbylevel);CHKERRQ(ierr);
71642eaa7f3SBarry Smith   for (j=1; j<=nlevels;j++) {
71742eaa7f3SBarry Smith     cnt = 0;
71842eaa7f3SBarry Smith     for (i=0; i<n; i++) {
71942eaa7f3SBarry Smith       if (mask[i]) continue;
72042eaa7f3SBarry Smith       if (level[i] != j) continue;
72142eaa7f3SBarry Smith       workid[cnt]         = i;
72242eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
72342eaa7f3SBarry Smith     }
72442eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
72542eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
72642eaa7f3SBarry Smith     ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr);
72742eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
72842eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
72942eaa7f3SBarry Smith     ierr = PetscMemcpy(idbylevel+tcnt,workid,cnt*sizeof(PetscInt));CHKERRQ(ierr);
73042eaa7f3SBarry Smith     tcnt += cnt;
73142eaa7f3SBarry Smith   }
73242eaa7f3SBarry Smith   if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes");
73342eaa7f3SBarry Smith   ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr);
73442eaa7f3SBarry Smith 
73542eaa7f3SBarry Smith   /* for each node list its column */
73642eaa7f3SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&column);CHKERRQ(ierr);
73742eaa7f3SBarry Smith   cnt = 0;
73842eaa7f3SBarry Smith   for (j=0; j<nlevels; j++) {
73942eaa7f3SBarry Smith     for (i=0; i<levelcnt[j]; i++) {
74042eaa7f3SBarry Smith       column[idbylevel[cnt++]] = i;
74142eaa7f3SBarry Smith     }
74242eaa7f3SBarry Smith   }
74342eaa7f3SBarry Smith 
74442eaa7f3SBarry Smith   *Nlevels   = nlevels;
74542eaa7f3SBarry Smith   *Level     = level;
74642eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
74742eaa7f3SBarry Smith   *Idbylevel = idbylevel;
74842eaa7f3SBarry Smith   *Column    = column;
74942eaa7f3SBarry Smith   PetscFunctionReturn(0);
75042eaa7f3SBarry Smith }
751