xref: /petsc/src/sys/utils/sorti.c (revision d7aa01a84ba4076bd64b3daf011b0beee9eee36c)
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 */
212*d7aa01a8SHong 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) {
220*d7aa01a8SHong 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   }
224*d7aa01a8SHong Zhang   SWAP3(L[0],L[right/2],J[0],J[right/2],K[0],K[right/2],tmp);
225*d7aa01a8SHong Zhang   vl   = L[0];
226c1f0200aSDmitry Karpeev   last = 0;
227c1f0200aSDmitry Karpeev   for (i=1; i<=right; i++) {
228*d7aa01a8SHong Zhang     if (L[i] < vl) {last++; SWAP3(L[last],L[i],J[last],J[i],K[last],K[i],tmp);}
229c1f0200aSDmitry Karpeev   }
230*d7aa01a8SHong Zhang   SWAP3(L[0],L[last],J[0],J[last],K[0],K[last],tmp);
231*d7aa01a8SHong Zhang   ierr = PetscSortIntWithArrayPair_Private(L,J,K,last-1);CHKERRQ(ierr);
232*d7aa01a8SHong 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 @*/
256*d7aa01a8SHong 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++) {
264*d7aa01a8SHong Zhang       ik = L[k];
265c1f0200aSDmitry Karpeev       for (j=k+1; j<n; j++) {
266*d7aa01a8SHong Zhang 	if (ik > L[j]) {
267*d7aa01a8SHong Zhang 	  SWAP3(L[k],L[j],J[k],J[j],K[k],K[j],tmp);
268*d7aa01a8SHong Zhang 	  ik = L[k];
269c1f0200aSDmitry Karpeev 	}
270c1f0200aSDmitry Karpeev       }
271c1f0200aSDmitry Karpeev     }
272c1f0200aSDmitry Karpeev   } else {
273*d7aa01a8SHong Zhang     ierr = PetscSortIntWithArrayPair_Private(L,J,K,n-1);CHKERRQ(ierr);
274c1f0200aSDmitry Karpeev   }
275c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
276c1f0200aSDmitry Karpeev }
277c1f0200aSDmitry Karpeev 
278c1f0200aSDmitry Karpeev 
2794d615ea0SBarry Smith #undef __FUNCT__
2804d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray_Private"
2814d615ea0SBarry Smith /*
2824d615ea0SBarry Smith    A simple version of quicksort; taken from Kernighan and Ritchie, page 87.
2834d615ea0SBarry Smith    Assumes 0 origin for v, number of elements = right+1 (right is index of
2844d615ea0SBarry Smith    right-most member).
2854d615ea0SBarry Smith */
2864d615ea0SBarry Smith static PetscErrorCode PetscSortMPIIntWithArray_Private(PetscMPIInt *v,PetscMPIInt *V,PetscMPIInt right)
2874d615ea0SBarry Smith {
2884d615ea0SBarry Smith   PetscErrorCode ierr;
2894d615ea0SBarry Smith   PetscMPIInt    i,vl,last,tmp;
2904d615ea0SBarry Smith 
2914d615ea0SBarry Smith   PetscFunctionBegin;
2924d615ea0SBarry Smith   if (right <= 1) {
2934d615ea0SBarry Smith     if (right == 1) {
2944d615ea0SBarry Smith       if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp);
2954d615ea0SBarry Smith     }
2964d615ea0SBarry Smith     PetscFunctionReturn(0);
2974d615ea0SBarry Smith   }
2984d615ea0SBarry Smith   SWAP2(v[0],v[right/2],V[0],V[right/2],tmp);
2994d615ea0SBarry Smith   vl   = v[0];
3004d615ea0SBarry Smith   last = 0;
3014d615ea0SBarry Smith   for (i=1; i<=right; i++) {
3024d615ea0SBarry Smith     if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);}
3034d615ea0SBarry Smith   }
3044d615ea0SBarry Smith   SWAP2(v[0],v[last],V[0],V[last],tmp);
3054d615ea0SBarry Smith   ierr = PetscSortMPIIntWithArray_Private(v,V,last-1);CHKERRQ(ierr);
3064d615ea0SBarry Smith   ierr = PetscSortMPIIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr);
3074d615ea0SBarry Smith   PetscFunctionReturn(0);
3084d615ea0SBarry Smith }
3094d615ea0SBarry Smith 
3104d615ea0SBarry Smith #undef __FUNCT__
3114d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray"
3124d615ea0SBarry Smith /*@
3134d615ea0SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order;
3144d615ea0SBarry Smith        changes a second array to match the sorted first array.
3154d615ea0SBarry Smith 
3164d615ea0SBarry Smith    Not Collective
3174d615ea0SBarry Smith 
3184d615ea0SBarry Smith    Input Parameters:
3194d615ea0SBarry Smith +  n  - number of values
3204d615ea0SBarry Smith .  i  - array of integers
3214d615ea0SBarry Smith -  I - second array of integers
3224d615ea0SBarry Smith 
3234d615ea0SBarry Smith    Level: intermediate
3244d615ea0SBarry Smith 
3254d615ea0SBarry Smith    Concepts: sorting^ints with array
3264d615ea0SBarry Smith 
3274d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
3284d615ea0SBarry Smith @*/
3297087cfbeSBarry Smith PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt i[],PetscMPIInt Ii[])
3304d615ea0SBarry Smith {
3314d615ea0SBarry Smith   PetscErrorCode ierr;
3324d615ea0SBarry Smith   PetscMPIInt    j,k,tmp,ik;
3334d615ea0SBarry Smith 
3344d615ea0SBarry Smith   PetscFunctionBegin;
3354d615ea0SBarry Smith   if (n<8) {
3364d615ea0SBarry Smith     for (k=0; k<n; k++) {
3374d615ea0SBarry Smith       ik = i[k];
3384d615ea0SBarry Smith       for (j=k+1; j<n; j++) {
3394d615ea0SBarry Smith 	if (ik > i[j]) {
3404d615ea0SBarry Smith 	  SWAP2(i[k],i[j],Ii[k],Ii[j],tmp);
3414d615ea0SBarry Smith 	  ik = i[k];
3424d615ea0SBarry Smith 	}
3434d615ea0SBarry Smith       }
3444d615ea0SBarry Smith     }
3454d615ea0SBarry Smith   } else {
3464d615ea0SBarry Smith     ierr = PetscSortMPIIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr);
3474d615ea0SBarry Smith   }
3484d615ea0SBarry Smith   PetscFunctionReturn(0);
3494d615ea0SBarry Smith }
3504d615ea0SBarry Smith 
351e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/
352e5c89e4eSSatish Balay #define SWAP2IntScalar(a,b,c,d,t,ts) {t=a;a=b;b=t;ts=c;c=d;d=ts;}
353e5c89e4eSSatish Balay 
354e5c89e4eSSatish Balay #undef __FUNCT__
355e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray_Private"
356e5c89e4eSSatish Balay /*
357e5c89e4eSSatish Balay    Modified from PetscSortIntWithArray_Private().
358e5c89e4eSSatish Balay */
359e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithScalarArray_Private(PetscInt *v,PetscScalar *V,PetscInt right)
360e5c89e4eSSatish Balay {
361e5c89e4eSSatish Balay   PetscErrorCode ierr;
362e5c89e4eSSatish Balay   PetscInt       i,vl,last,tmp;
363e5c89e4eSSatish Balay   PetscScalar    stmp;
364e5c89e4eSSatish Balay 
365e5c89e4eSSatish Balay   PetscFunctionBegin;
366e5c89e4eSSatish Balay   if (right <= 1) {
367e5c89e4eSSatish Balay     if (right == 1) {
368e5c89e4eSSatish Balay       if (v[0] > v[1]) SWAP2IntScalar(v[0],v[1],V[0],V[1],tmp,stmp);
369e5c89e4eSSatish Balay     }
370e5c89e4eSSatish Balay     PetscFunctionReturn(0);
371e5c89e4eSSatish Balay   }
372e5c89e4eSSatish Balay   SWAP2IntScalar(v[0],v[right/2],V[0],V[right/2],tmp,stmp);
373e5c89e4eSSatish Balay   vl   = v[0];
374e5c89e4eSSatish Balay   last = 0;
375e5c89e4eSSatish Balay   for (i=1; i<=right; i++) {
376e5c89e4eSSatish Balay     if (v[i] < vl) {last++; SWAP2IntScalar(v[last],v[i],V[last],V[i],tmp,stmp);}
377e5c89e4eSSatish Balay   }
378e5c89e4eSSatish Balay   SWAP2IntScalar(v[0],v[last],V[0],V[last],tmp,stmp);
379e5c89e4eSSatish Balay   ierr = PetscSortIntWithScalarArray_Private(v,V,last-1);CHKERRQ(ierr);
380e5c89e4eSSatish Balay   ierr = PetscSortIntWithScalarArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr);
381e5c89e4eSSatish Balay   PetscFunctionReturn(0);
382e5c89e4eSSatish Balay }
383e5c89e4eSSatish Balay 
384e5c89e4eSSatish Balay #undef __FUNCT__
385e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray"
386e5c89e4eSSatish Balay /*@
387e5c89e4eSSatish Balay    PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
388e5c89e4eSSatish Balay        changes a second SCALAR array to match the sorted first INTEGER array.
389e5c89e4eSSatish Balay 
390e5c89e4eSSatish Balay    Not Collective
391e5c89e4eSSatish Balay 
392e5c89e4eSSatish Balay    Input Parameters:
393e5c89e4eSSatish Balay +  n  - number of values
394e5c89e4eSSatish Balay .  i  - array of integers
395e5c89e4eSSatish Balay -  I - second array of scalars
396e5c89e4eSSatish Balay 
397e5c89e4eSSatish Balay    Level: intermediate
398e5c89e4eSSatish Balay 
399e5c89e4eSSatish Balay    Concepts: sorting^ints with array
400e5c89e4eSSatish Balay 
401e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
402e5c89e4eSSatish Balay @*/
4037087cfbeSBarry Smith PetscErrorCode  PetscSortIntWithScalarArray(PetscInt n,PetscInt i[],PetscScalar Ii[])
404e5c89e4eSSatish Balay {
405e5c89e4eSSatish Balay   PetscErrorCode ierr;
406e5c89e4eSSatish Balay   PetscInt       j,k,tmp,ik;
407e5c89e4eSSatish Balay   PetscScalar    stmp;
408e5c89e4eSSatish Balay 
409e5c89e4eSSatish Balay   PetscFunctionBegin;
410e5c89e4eSSatish Balay   if (n<8) {
411e5c89e4eSSatish Balay     for (k=0; k<n; k++) {
412e5c89e4eSSatish Balay       ik = i[k];
413e5c89e4eSSatish Balay       for (j=k+1; j<n; j++) {
414e5c89e4eSSatish Balay 	if (ik > i[j]) {
415b7940d39SSatish Balay 	  SWAP2IntScalar(i[k],i[j],Ii[k],Ii[j],tmp,stmp);
416e5c89e4eSSatish Balay 	  ik = i[k];
417e5c89e4eSSatish Balay 	}
418e5c89e4eSSatish Balay       }
419e5c89e4eSSatish Balay     }
420e5c89e4eSSatish Balay   } else {
421b7940d39SSatish Balay     ierr = PetscSortIntWithScalarArray_Private(i,Ii,n-1);CHKERRQ(ierr);
422e5c89e4eSSatish Balay   }
423e5c89e4eSSatish Balay   PetscFunctionReturn(0);
424e5c89e4eSSatish Balay }
425e5c89e4eSSatish Balay 
426c1f0200aSDmitry Karpeev #undef __FUNCT__
427c1f0200aSDmitry Karpeev #define __FUNCT__ "PetscMergeIntArrayPair"
428c1f0200aSDmitry Karpeev /*@
429c1f0200aSDmitry Karpeev    PetscMergeIntArrayPair -     Merges two SORTED integer arrays along with an additional array of integers.
430c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
431c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
432c1f0200aSDmitry Karpeev 
433c1f0200aSDmitry Karpeev    Not Collective
434c1f0200aSDmitry Karpeev 
435c1f0200aSDmitry Karpeev    Input Parameters:
436c1f0200aSDmitry Karpeev +  an  - number of values in the first array
437c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
438c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
439c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
440c1f0200aSDmitry Karpeev .  bI  - second array of integers
441c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
442c1f0200aSDmitry Karpeev 
443c1f0200aSDmitry Karpeev    Output Parameters:
444c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
445c1f0200aSDmitry Karpeev .  I   - merged sorted array
446c1f0200aSDmitry Karpeev -  J   - merged additional array
447c1f0200aSDmitry Karpeev 
448c1f0200aSDmitry Karpeev    Level: intermediate
449c1f0200aSDmitry Karpeev 
450c1f0200aSDmitry Karpeev    Concepts: merging^arrays
451c1f0200aSDmitry Karpeev 
452c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
453c1f0200aSDmitry Karpeev @*/
454*d7aa01a8SHong 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)
455c1f0200aSDmitry Karpeev {
456c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
457*d7aa01a8SHong Zhang   PetscInt n_, *L_ = *L, *J_= *J, ak, bk, k;
458c1f0200aSDmitry Karpeev 
459c1f0200aSDmitry Karpeev   n_ = an + bn;
460c1f0200aSDmitry Karpeev   *n = n_;
461*d7aa01a8SHong Zhang   if(!L_) {
462*d7aa01a8SHong Zhang     ierr = PetscMalloc(n_*sizeof(PetscInt), L); CHKERRQ(ierr);
463*d7aa01a8SHong Zhang     L_ = *L;
464c1f0200aSDmitry Karpeev   }
465c1f0200aSDmitry Karpeev   if(!J_){
466c1f0200aSDmitry Karpeev     ierr = PetscMalloc(n_*sizeof(PetscInt), &J_); CHKERRQ(ierr);
467c1f0200aSDmitry Karpeev     J_ = *J;
468c1f0200aSDmitry Karpeev   }
469c1f0200aSDmitry Karpeev   k = ak = bk = 0;
470c1f0200aSDmitry Karpeev   while(ak < an && bk < bn) {
471c1f0200aSDmitry Karpeev     if(aI[ak] <= bI[bk]) {
472*d7aa01a8SHong Zhang       L_[k] = aI[ak];
473c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
474c1f0200aSDmitry Karpeev       ++ak;
475c1f0200aSDmitry Karpeev       ++k;
476c1f0200aSDmitry Karpeev     }
477c1f0200aSDmitry Karpeev     else {
478*d7aa01a8SHong Zhang       L_[k] = bI[bk];
479c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
480c1f0200aSDmitry Karpeev       ++bk;
481c1f0200aSDmitry Karpeev       ++k;
482c1f0200aSDmitry Karpeev     }
483c1f0200aSDmitry Karpeev   }
484c1f0200aSDmitry Karpeev   if(ak < an) {
485*d7aa01a8SHong Zhang     ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt)); CHKERRQ(ierr);
486c1f0200aSDmitry Karpeev     ierr = PetscMemcpy(J_+k,aJ+ak,(an-ak)*sizeof(PetscInt)); CHKERRQ(ierr);
487c1f0200aSDmitry Karpeev     k += (an-ak);
488c1f0200aSDmitry Karpeev   }
489c1f0200aSDmitry Karpeev   if(bk < bn) {
490*d7aa01a8SHong Zhang     ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt)); CHKERRQ(ierr);
491c1f0200aSDmitry Karpeev     ierr = PetscMemcpy(J_+k,bJ+bk,(bn-bk)*sizeof(PetscInt)); CHKERRQ(ierr);
492c1f0200aSDmitry Karpeev     k += (bn-bk);
493c1f0200aSDmitry Karpeev   }
494c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
495c1f0200aSDmitry Karpeev }
496c1f0200aSDmitry Karpeev 
497e5c89e4eSSatish Balay 
49842eaa7f3SBarry Smith #undef __FUNCT__
49942eaa7f3SBarry Smith #define __FUNCT__ "PetscProcessTree"
50042eaa7f3SBarry Smith /*@
50142eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
50242eaa7f3SBarry Smith 
50342eaa7f3SBarry Smith    Not Collective
50442eaa7f3SBarry Smith 
50542eaa7f3SBarry Smith    Input Parameters:
50642eaa7f3SBarry Smith +  n  - number of values
50742eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
50842eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
50942eaa7f3SBarry Smith 
51042eaa7f3SBarry Smith    Output Parameters:
51142eaa7f3SBarry Smith +  Nlevels - the number of levels
51242eaa7f3SBarry Smith .  Level - for each node tells its level
51342eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
51442eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
51542eaa7f3SBarry Smith -  Column - for each id tells its column index
51642eaa7f3SBarry Smith 
51742eaa7f3SBarry Smith    Level: intermediate
51842eaa7f3SBarry Smith 
51942eaa7f3SBarry Smith 
52042eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation()
52142eaa7f3SBarry Smith @*/
5227087cfbeSBarry Smith PetscErrorCode  PetscProcessTree(PetscInt n,const PetscBool  mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column)
52342eaa7f3SBarry Smith {
52442eaa7f3SBarry Smith   PetscInt       i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column;
52542eaa7f3SBarry Smith   PetscErrorCode ierr;
526ace3abfcSBarry Smith   PetscBool      done = PETSC_FALSE;
52742eaa7f3SBarry Smith 
52842eaa7f3SBarry Smith   PetscFunctionBegin;
52942eaa7f3SBarry Smith   if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set");
53042eaa7f3SBarry Smith   for (i=0; i<n; i++) {
53142eaa7f3SBarry Smith     if (mask[i]) continue;
53242eaa7f3SBarry Smith     if (parentid[i]  == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent");
53342eaa7f3SBarry Smith     if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked");
53442eaa7f3SBarry Smith   }
53542eaa7f3SBarry Smith 
53642eaa7f3SBarry Smith 
53742eaa7f3SBarry Smith   for (i=0; i<n; i++) {
53842eaa7f3SBarry Smith     if (!mask[i]) nmask++;
53942eaa7f3SBarry Smith   }
54042eaa7f3SBarry Smith 
54142eaa7f3SBarry Smith   /* determine the level in the tree of each node */
54242eaa7f3SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&level);CHKERRQ(ierr);
54342eaa7f3SBarry Smith   ierr = PetscMemzero(level,n*sizeof(PetscInt));CHKERRQ(ierr);
54442eaa7f3SBarry Smith   level[0] = 1;
54542eaa7f3SBarry Smith   while (!done) {
54642eaa7f3SBarry Smith     done = PETSC_TRUE;
54742eaa7f3SBarry Smith     for (i=0; i<n; i++) {
54842eaa7f3SBarry Smith       if (mask[i]) continue;
54942eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
55042eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
55142eaa7f3SBarry Smith     }
55242eaa7f3SBarry Smith   }
55342eaa7f3SBarry Smith   for (i=0; i<n; i++) {
55442eaa7f3SBarry Smith     level[i]--;
55542eaa7f3SBarry Smith     nlevels = PetscMax(nlevels,level[i]);
55642eaa7f3SBarry Smith   }
55742eaa7f3SBarry Smith 
55842eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
55942eaa7f3SBarry Smith   ierr = PetscMalloc(nlevels*sizeof(PetscInt),&levelcnt);CHKERRQ(ierr);
56042eaa7f3SBarry Smith   ierr = PetscMemzero(levelcnt,nlevels*sizeof(PetscInt));CHKERRQ(ierr);
56142eaa7f3SBarry Smith   for (i=0; i<n; i++) {
56242eaa7f3SBarry Smith     if (mask[i]) continue;
56342eaa7f3SBarry Smith     levelcnt[level[i]-1]++;
56442eaa7f3SBarry Smith   }
56542eaa7f3SBarry Smith   for (i=0; i<nlevels;i++) {
56642eaa7f3SBarry Smith     levelmax = PetscMax(levelmax,levelcnt[i]);
56742eaa7f3SBarry Smith   }
56842eaa7f3SBarry Smith 
56942eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
57042eaa7f3SBarry Smith   ierr = PetscMalloc2(levelmax,PetscInt,&workid,levelmax,PetscInt,&workparentid);CHKERRQ(ierr);
57142eaa7f3SBarry Smith   ierr = PetscMalloc(nmask*sizeof(PetscInt),&idbylevel);CHKERRQ(ierr);
57242eaa7f3SBarry Smith   for (j=1; j<=nlevels;j++) {
57342eaa7f3SBarry Smith     cnt = 0;
57442eaa7f3SBarry Smith     for (i=0; i<n; i++) {
57542eaa7f3SBarry Smith       if (mask[i]) continue;
57642eaa7f3SBarry Smith       if (level[i] != j) continue;
57742eaa7f3SBarry Smith       workid[cnt]         = i;
57842eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
57942eaa7f3SBarry Smith     }
58042eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
58142eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
58242eaa7f3SBarry Smith     ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr);
58342eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
58442eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
58542eaa7f3SBarry Smith     ierr = PetscMemcpy(idbylevel+tcnt,workid,cnt*sizeof(PetscInt));CHKERRQ(ierr);
58642eaa7f3SBarry Smith     tcnt += cnt;
58742eaa7f3SBarry Smith   }
58842eaa7f3SBarry Smith   if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes");
58942eaa7f3SBarry Smith   ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr);
59042eaa7f3SBarry Smith 
59142eaa7f3SBarry Smith   /* for each node list its column */
59242eaa7f3SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&column);CHKERRQ(ierr);
59342eaa7f3SBarry Smith   cnt = 0;
59442eaa7f3SBarry Smith   for (j=0; j<nlevels; j++) {
59542eaa7f3SBarry Smith     for (i=0; i<levelcnt[j]; i++) {
59642eaa7f3SBarry Smith       column[idbylevel[cnt++]] = i;
59742eaa7f3SBarry Smith     }
59842eaa7f3SBarry Smith   }
59942eaa7f3SBarry Smith 
60042eaa7f3SBarry Smith   *Nlevels   = nlevels;
60142eaa7f3SBarry Smith   *Level     = level;
60242eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
60342eaa7f3SBarry Smith   *Idbylevel = idbylevel;
60442eaa7f3SBarry Smith   *Column    = column;
60542eaa7f3SBarry Smith   PetscFunctionReturn(0);
60642eaa7f3SBarry Smith }
607