xref: /petsc/src/sys/utils/sorti.c (revision 17df18f2f5feaf86b818dab91272a58888dbc575)
17d0a6c19SBarry Smith 
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file contains routines for sorting integers. Values are sorted in place.
4e5c89e4eSSatish Balay  */
5afcb2eb5SJed Brown #include <petsc-private/petscimpl.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     }
85a297a907SKarl Rupp   } else PetscSortInt_Private(i,n-1);
86e5c89e4eSSatish Balay   PetscFunctionReturn(0);
87e5c89e4eSSatish Balay }
88e5c89e4eSSatish Balay 
891db36a52SBarry Smith #undef __FUNCT__
901db36a52SBarry Smith #define __FUNCT__ "PetscSortRemoveDupsInt"
911db36a52SBarry Smith /*@
921db36a52SBarry Smith    PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries
931db36a52SBarry Smith 
941db36a52SBarry Smith    Not Collective
951db36a52SBarry Smith 
961db36a52SBarry Smith    Input Parameters:
971db36a52SBarry Smith +  n  - number of values
981db36a52SBarry Smith -  ii  - array of integers
991db36a52SBarry Smith 
1001db36a52SBarry Smith    Output Parameter:
1011db36a52SBarry Smith .  n - number of non-redundant values
1021db36a52SBarry Smith 
1031db36a52SBarry Smith    Level: intermediate
1041db36a52SBarry Smith 
1051db36a52SBarry Smith    Concepts: sorting^ints
1061db36a52SBarry Smith 
1071db36a52SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
1081db36a52SBarry Smith @*/
1097087cfbeSBarry Smith PetscErrorCode  PetscSortRemoveDupsInt(PetscInt *n,PetscInt ii[])
1101db36a52SBarry Smith {
1111db36a52SBarry Smith   PetscErrorCode ierr;
1121db36a52SBarry Smith   PetscInt       i,s = 0,N = *n, b = 0;
1131db36a52SBarry Smith 
1141db36a52SBarry Smith   PetscFunctionBegin;
1151db36a52SBarry Smith   ierr = PetscSortInt(N,ii);CHKERRQ(ierr);
1161db36a52SBarry Smith   for (i=0; i<N-1; i++) {
117a297a907SKarl Rupp     if (ii[b+s+1] != ii[b]) {
118a297a907SKarl Rupp       ii[b+1] = ii[b+s+1]; b++;
119a297a907SKarl Rupp     } else s++;
1201db36a52SBarry Smith   }
1211db36a52SBarry Smith   *n = N - s;
1221db36a52SBarry Smith   PetscFunctionReturn(0);
1231db36a52SBarry Smith }
1241db36a52SBarry Smith 
12560e03357SMatthew G Knepley #undef __FUNCT__
12660e03357SMatthew G Knepley #define __FUNCT__ "PetscFindInt"
12760e03357SMatthew G Knepley /*@
128d9f1162dSJed Brown   PetscFindInt - Finds integer in a sorted array of integers
12960e03357SMatthew G Knepley 
13060e03357SMatthew G Knepley    Not Collective
13160e03357SMatthew G Knepley 
13260e03357SMatthew G Knepley    Input Parameters:
13360e03357SMatthew G Knepley +  key - the integer to locate
134d9f1162dSJed Brown .  n   - number of values in the array
13560e03357SMatthew G Knepley -  ii  - array of integers
13660e03357SMatthew G Knepley 
13760e03357SMatthew G Knepley    Output Parameter:
138d9f1162dSJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
13960e03357SMatthew G Knepley 
14060e03357SMatthew G Knepley    Level: intermediate
14160e03357SMatthew G Knepley 
14260e03357SMatthew G Knepley    Concepts: sorting^ints
14360e03357SMatthew G Knepley 
14460e03357SMatthew G Knepley .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
14560e03357SMatthew G Knepley @*/
146026ec6cbSMatthew G Knepley PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt ii[], PetscInt *loc)
14760e03357SMatthew G Knepley {
148d9f1162dSJed Brown   PetscInt lo = 0,hi = n;
14960e03357SMatthew G Knepley 
15060e03357SMatthew G Knepley   PetscFunctionBegin;
15160e03357SMatthew G Knepley   PetscValidPointer(loc,4);
15298781d82SMatthew G Knepley   if (!n) {*loc = -1; PetscFunctionReturn(0);}
15398781d82SMatthew G Knepley   PetscValidPointer(ii,3);
154d9f1162dSJed Brown   while (hi - lo > 1) {
155d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo)/2;
156d9f1162dSJed Brown     if (key < ii[mid]) hi = mid;
157d9f1162dSJed Brown     else               lo = mid;
15860e03357SMatthew G Knepley   }
159d9f1162dSJed 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     }
380a297a907SKarl Rupp   } else PetscSortMPIInt_Private(i,n-1);
38117d7d925SStefano Zampini   PetscFunctionReturn(0);
38217d7d925SStefano Zampini }
38317d7d925SStefano Zampini 
38417d7d925SStefano Zampini #undef __FUNCT__
38517d7d925SStefano Zampini #define __FUNCT__ "PetscSortRemoveDupsMPIInt"
38617d7d925SStefano Zampini /*@
38717d7d925SStefano Zampini    PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries
38817d7d925SStefano Zampini 
38917d7d925SStefano Zampini    Not Collective
39017d7d925SStefano Zampini 
39117d7d925SStefano Zampini    Input Parameters:
39217d7d925SStefano Zampini +  n  - number of values
39317d7d925SStefano Zampini -  ii  - array of integers
39417d7d925SStefano Zampini 
39517d7d925SStefano Zampini    Output Parameter:
39617d7d925SStefano Zampini .  n - number of non-redundant values
39717d7d925SStefano Zampini 
39817d7d925SStefano Zampini    Level: intermediate
39917d7d925SStefano Zampini 
40017d7d925SStefano Zampini    Concepts: sorting^ints
40117d7d925SStefano Zampini 
40217d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
40317d7d925SStefano Zampini @*/
40417d7d925SStefano Zampini PetscErrorCode  PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt ii[])
40517d7d925SStefano Zampini {
40617d7d925SStefano Zampini   PetscErrorCode ierr;
40717d7d925SStefano Zampini   PetscInt       i,s = 0,N = *n, b = 0;
40817d7d925SStefano Zampini 
40917d7d925SStefano Zampini   PetscFunctionBegin;
41017d7d925SStefano Zampini   ierr = PetscSortMPIInt(N,ii);CHKERRQ(ierr);
41117d7d925SStefano Zampini   for (i=0; i<N-1; i++) {
412a297a907SKarl Rupp     if (ii[b+s+1] != ii[b]) {
413a297a907SKarl Rupp       ii[b+1] = ii[b+s+1]; b++;
414a297a907SKarl Rupp     } else s++;
41517d7d925SStefano Zampini   }
41617d7d925SStefano Zampini   *n = N - s;
41717d7d925SStefano Zampini   PetscFunctionReturn(0);
41817d7d925SStefano Zampini }
419c1f0200aSDmitry Karpeev 
4204d615ea0SBarry Smith #undef __FUNCT__
4214d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray_Private"
4224d615ea0SBarry Smith /*
4234d615ea0SBarry Smith    A simple version of quicksort; taken from Kernighan and Ritchie, page 87.
4244d615ea0SBarry Smith    Assumes 0 origin for v, number of elements = right+1 (right is index of
4254d615ea0SBarry Smith    right-most member).
4264d615ea0SBarry Smith */
4274d615ea0SBarry Smith static PetscErrorCode PetscSortMPIIntWithArray_Private(PetscMPIInt *v,PetscMPIInt *V,PetscMPIInt right)
4284d615ea0SBarry Smith {
4294d615ea0SBarry Smith   PetscErrorCode ierr;
4304d615ea0SBarry Smith   PetscMPIInt    i,vl,last,tmp;
4314d615ea0SBarry Smith 
4324d615ea0SBarry Smith   PetscFunctionBegin;
4334d615ea0SBarry Smith   if (right <= 1) {
4344d615ea0SBarry Smith     if (right == 1) {
4354d615ea0SBarry Smith       if (v[0] > v[1]) SWAP2(v[0],v[1],V[0],V[1],tmp);
4364d615ea0SBarry Smith     }
4374d615ea0SBarry Smith     PetscFunctionReturn(0);
4384d615ea0SBarry Smith   }
4394d615ea0SBarry Smith   SWAP2(v[0],v[right/2],V[0],V[right/2],tmp);
4404d615ea0SBarry Smith   vl   = v[0];
4414d615ea0SBarry Smith   last = 0;
4424d615ea0SBarry Smith   for (i=1; i<=right; i++) {
4434d615ea0SBarry Smith     if (v[i] < vl) {last++; SWAP2(v[last],v[i],V[last],V[i],tmp);}
4444d615ea0SBarry Smith   }
4454d615ea0SBarry Smith   SWAP2(v[0],v[last],V[0],V[last],tmp);
4464d615ea0SBarry Smith   ierr = PetscSortMPIIntWithArray_Private(v,V,last-1);CHKERRQ(ierr);
4474d615ea0SBarry Smith   ierr = PetscSortMPIIntWithArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr);
4484d615ea0SBarry Smith   PetscFunctionReturn(0);
4494d615ea0SBarry Smith }
4504d615ea0SBarry Smith 
4514d615ea0SBarry Smith #undef __FUNCT__
4524d615ea0SBarry Smith #define __FUNCT__ "PetscSortMPIIntWithArray"
4534d615ea0SBarry Smith /*@
4544d615ea0SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order;
4554d615ea0SBarry Smith        changes a second array to match the sorted first array.
4564d615ea0SBarry Smith 
4574d615ea0SBarry Smith    Not Collective
4584d615ea0SBarry Smith 
4594d615ea0SBarry Smith    Input Parameters:
4604d615ea0SBarry Smith +  n  - number of values
4614d615ea0SBarry Smith .  i  - array of integers
4624d615ea0SBarry Smith -  I - second array of integers
4634d615ea0SBarry Smith 
4644d615ea0SBarry Smith    Level: intermediate
4654d615ea0SBarry Smith 
4664d615ea0SBarry Smith    Concepts: sorting^ints with array
4674d615ea0SBarry Smith 
4684d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
4694d615ea0SBarry Smith @*/
4707087cfbeSBarry Smith PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt i[],PetscMPIInt Ii[])
4714d615ea0SBarry Smith {
4724d615ea0SBarry Smith   PetscErrorCode ierr;
4734d615ea0SBarry Smith   PetscMPIInt    j,k,tmp,ik;
4744d615ea0SBarry Smith 
4754d615ea0SBarry Smith   PetscFunctionBegin;
4764d615ea0SBarry Smith   if (n<8) {
4774d615ea0SBarry Smith     for (k=0; k<n; k++) {
4784d615ea0SBarry Smith       ik = i[k];
4794d615ea0SBarry Smith       for (j=k+1; j<n; j++) {
4804d615ea0SBarry Smith         if (ik > i[j]) {
4814d615ea0SBarry Smith           SWAP2(i[k],i[j],Ii[k],Ii[j],tmp);
4824d615ea0SBarry Smith           ik = i[k];
4834d615ea0SBarry Smith         }
4844d615ea0SBarry Smith       }
4854d615ea0SBarry Smith     }
4864d615ea0SBarry Smith   } else {
4874d615ea0SBarry Smith     ierr = PetscSortMPIIntWithArray_Private(i,Ii,n-1);CHKERRQ(ierr);
4884d615ea0SBarry Smith   }
4894d615ea0SBarry Smith   PetscFunctionReturn(0);
4904d615ea0SBarry Smith }
4914d615ea0SBarry Smith 
492e5c89e4eSSatish Balay /* -----------------------------------------------------------------------*/
493e5c89e4eSSatish Balay #define SWAP2IntScalar(a,b,c,d,t,ts) {t=a;a=b;b=t;ts=c;c=d;d=ts;}
494e5c89e4eSSatish Balay 
495e5c89e4eSSatish Balay #undef __FUNCT__
496e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray_Private"
497e5c89e4eSSatish Balay /*
498e5c89e4eSSatish Balay    Modified from PetscSortIntWithArray_Private().
499e5c89e4eSSatish Balay */
500e5c89e4eSSatish Balay static PetscErrorCode PetscSortIntWithScalarArray_Private(PetscInt *v,PetscScalar *V,PetscInt right)
501e5c89e4eSSatish Balay {
502e5c89e4eSSatish Balay   PetscErrorCode ierr;
503e5c89e4eSSatish Balay   PetscInt       i,vl,last,tmp;
504e5c89e4eSSatish Balay   PetscScalar    stmp;
505e5c89e4eSSatish Balay 
506e5c89e4eSSatish Balay   PetscFunctionBegin;
507e5c89e4eSSatish Balay   if (right <= 1) {
508e5c89e4eSSatish Balay     if (right == 1) {
509e5c89e4eSSatish Balay       if (v[0] > v[1]) SWAP2IntScalar(v[0],v[1],V[0],V[1],tmp,stmp);
510e5c89e4eSSatish Balay     }
511e5c89e4eSSatish Balay     PetscFunctionReturn(0);
512e5c89e4eSSatish Balay   }
513e5c89e4eSSatish Balay   SWAP2IntScalar(v[0],v[right/2],V[0],V[right/2],tmp,stmp);
514e5c89e4eSSatish Balay   vl   = v[0];
515e5c89e4eSSatish Balay   last = 0;
516e5c89e4eSSatish Balay   for (i=1; i<=right; i++) {
517e5c89e4eSSatish Balay     if (v[i] < vl) {last++; SWAP2IntScalar(v[last],v[i],V[last],V[i],tmp,stmp);}
518e5c89e4eSSatish Balay   }
519e5c89e4eSSatish Balay   SWAP2IntScalar(v[0],v[last],V[0],V[last],tmp,stmp);
520e5c89e4eSSatish Balay   ierr = PetscSortIntWithScalarArray_Private(v,V,last-1);CHKERRQ(ierr);
521e5c89e4eSSatish Balay   ierr = PetscSortIntWithScalarArray_Private(v+last+1,V+last+1,right-(last+1));CHKERRQ(ierr);
522e5c89e4eSSatish Balay   PetscFunctionReturn(0);
523e5c89e4eSSatish Balay }
524e5c89e4eSSatish Balay 
525e5c89e4eSSatish Balay #undef __FUNCT__
526e5c89e4eSSatish Balay #define __FUNCT__ "PetscSortIntWithScalarArray"
527e5c89e4eSSatish Balay /*@
528e5c89e4eSSatish Balay    PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
529e5c89e4eSSatish Balay        changes a second SCALAR array to match the sorted first INTEGER array.
530e5c89e4eSSatish Balay 
531e5c89e4eSSatish Balay    Not Collective
532e5c89e4eSSatish Balay 
533e5c89e4eSSatish Balay    Input Parameters:
534e5c89e4eSSatish Balay +  n  - number of values
535e5c89e4eSSatish Balay .  i  - array of integers
536e5c89e4eSSatish Balay -  I - second array of scalars
537e5c89e4eSSatish Balay 
538e5c89e4eSSatish Balay    Level: intermediate
539e5c89e4eSSatish Balay 
540e5c89e4eSSatish Balay    Concepts: sorting^ints with array
541e5c89e4eSSatish Balay 
542e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
543e5c89e4eSSatish Balay @*/
5447087cfbeSBarry Smith PetscErrorCode  PetscSortIntWithScalarArray(PetscInt n,PetscInt i[],PetscScalar Ii[])
545e5c89e4eSSatish Balay {
546e5c89e4eSSatish Balay   PetscErrorCode ierr;
547e5c89e4eSSatish Balay   PetscInt       j,k,tmp,ik;
548e5c89e4eSSatish Balay   PetscScalar    stmp;
549e5c89e4eSSatish Balay 
550e5c89e4eSSatish Balay   PetscFunctionBegin;
551e5c89e4eSSatish Balay   if (n<8) {
552e5c89e4eSSatish Balay     for (k=0; k<n; k++) {
553e5c89e4eSSatish Balay       ik = i[k];
554e5c89e4eSSatish Balay       for (j=k+1; j<n; j++) {
555e5c89e4eSSatish Balay         if (ik > i[j]) {
556b7940d39SSatish Balay           SWAP2IntScalar(i[k],i[j],Ii[k],Ii[j],tmp,stmp);
557e5c89e4eSSatish Balay           ik = i[k];
558e5c89e4eSSatish Balay         }
559e5c89e4eSSatish Balay       }
560e5c89e4eSSatish Balay     }
561e5c89e4eSSatish Balay   } else {
562b7940d39SSatish Balay     ierr = PetscSortIntWithScalarArray_Private(i,Ii,n-1);CHKERRQ(ierr);
563e5c89e4eSSatish Balay   }
564e5c89e4eSSatish Balay   PetscFunctionReturn(0);
565e5c89e4eSSatish Balay }
566e5c89e4eSSatish Balay 
567*17df18f2SToby Isaac #define SWAP2IntData(a,b,c,d,t,td,siz) {t=a;a=b;b=t;memcpy(td,c,siz);memcpy(c,d,siz);memcpy(d,td,siz);}
568*17df18f2SToby Isaac 
569*17df18f2SToby Isaac #undef __FUNCT__
570*17df18f2SToby Isaac #define __FUNCT__ "PetscSortIntWithDataArray_Private"
571*17df18f2SToby Isaac /*
572*17df18f2SToby Isaac    Modified from PetscSortIntWithArray_Private().
573*17df18f2SToby Isaac */
574*17df18f2SToby Isaac static PetscErrorCode PetscSortIntWithDataArray_Private(PetscInt *v,char *V,PetscInt right,size_t size,void *work)
575*17df18f2SToby Isaac {
576*17df18f2SToby Isaac   PetscErrorCode ierr;
577*17df18f2SToby Isaac   PetscInt       i,vl,last,tmp;
578*17df18f2SToby Isaac 
579*17df18f2SToby Isaac   PetscFunctionBegin;
580*17df18f2SToby Isaac   if (right <= 1) {
581*17df18f2SToby Isaac     if (right == 1) {
582*17df18f2SToby Isaac       if (v[0] > v[1]) SWAP2IntData(v[0],v[1],V,V+size,tmp,work,size);
583*17df18f2SToby Isaac     }
584*17df18f2SToby Isaac     PetscFunctionReturn(0);
585*17df18f2SToby Isaac   }
586*17df18f2SToby Isaac   SWAP2IntData(v[0],v[right/2],V,V+size*(right/2),tmp,work,size);
587*17df18f2SToby Isaac   vl   = v[0];
588*17df18f2SToby Isaac   last = 0;
589*17df18f2SToby Isaac   for (i=1; i<=right; i++) {
590*17df18f2SToby Isaac     if (v[i] < vl) {last++; SWAP2IntData(v[last],v[i],V+size*last,V+size*i,tmp,work,size);}
591*17df18f2SToby Isaac   }
592*17df18f2SToby Isaac   SWAP2IntData(v[0],v[last],V,V + size*last,tmp,work,size);
593*17df18f2SToby Isaac   ierr = PetscSortIntWithDataArray_Private(v,V,last-1,size,work);CHKERRQ(ierr);
594*17df18f2SToby Isaac   ierr = PetscSortIntWithDataArray_Private(v+last+1,V+size*(last+1),right-(last+1),size,work);CHKERRQ(ierr);
595*17df18f2SToby Isaac   PetscFunctionReturn(0);
596*17df18f2SToby Isaac }
597*17df18f2SToby Isaac 
598*17df18f2SToby Isaac #undef __FUNCT__
599*17df18f2SToby Isaac #define __FUNCT__ "PetscSortIntWithDataArray"
600*17df18f2SToby Isaac /*@
601*17df18f2SToby Isaac    PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order;
602*17df18f2SToby Isaac        changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
603*17df18f2SToby Isaac        provide workspace (the size of an element in the data array) to use when sorting.
604*17df18f2SToby Isaac 
605*17df18f2SToby Isaac    Not Collective
606*17df18f2SToby Isaac 
607*17df18f2SToby Isaac    Input Parameters:
608*17df18f2SToby Isaac +  n  - number of values
609*17df18f2SToby Isaac .  i  - array of integers
610*17df18f2SToby Isaac .  Ii - second array of data
611*17df18f2SToby Isaac .  size - sizeof elements in the data array in bytes
612*17df18f2SToby Isaac -  work - workspace of "size" bytes used when sorting
613*17df18f2SToby Isaac 
614*17df18f2SToby Isaac    Level: intermediate
615*17df18f2SToby Isaac 
616*17df18f2SToby Isaac    Concepts: sorting^ints with array
617*17df18f2SToby Isaac 
618*17df18f2SToby Isaac .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
619*17df18f2SToby Isaac @*/
620*17df18f2SToby Isaac PetscErrorCode  PetscSortIntWithDataArray(PetscInt n,PetscInt i[],void *Ii,size_t size,void *work)
621*17df18f2SToby Isaac {
622*17df18f2SToby Isaac   char           *V = (char *) Ii;
623*17df18f2SToby Isaac   PetscErrorCode ierr;
624*17df18f2SToby Isaac   PetscInt       j,k,tmp,ik;
625*17df18f2SToby Isaac 
626*17df18f2SToby Isaac   PetscFunctionBegin;
627*17df18f2SToby Isaac   if (n<8) {
628*17df18f2SToby Isaac     for (k=0; k<n; k++) {
629*17df18f2SToby Isaac       ik = i[k];
630*17df18f2SToby Isaac       for (j=k+1; j<n; j++) {
631*17df18f2SToby Isaac         if (ik > i[j]) {
632*17df18f2SToby Isaac           SWAP2IntData(i[k],i[j],V+size*k,V+size*j,tmp,work,size);
633*17df18f2SToby Isaac           ik = i[k];
634*17df18f2SToby Isaac         }
635*17df18f2SToby Isaac       }
636*17df18f2SToby Isaac     }
637*17df18f2SToby Isaac   } else {
638*17df18f2SToby Isaac     ierr = PetscSortIntWithDataArray_Private(i,V,n-1,size,work);CHKERRQ(ierr);
639*17df18f2SToby Isaac   }
640*17df18f2SToby Isaac   PetscFunctionReturn(0);
641*17df18f2SToby Isaac }
642*17df18f2SToby Isaac 
643b4301105SBarry Smith 
64421e72a00SBarry Smith #undef __FUNCT__
64521e72a00SBarry Smith #define __FUNCT__ "PetscMergeIntArray"
64621e72a00SBarry Smith /*@
64721e72a00SBarry Smith    PetscMergeIntArray -     Merges two SORTED integer arrays, removes duplicate elements.
64821e72a00SBarry Smith 
64921e72a00SBarry Smith    Not Collective
65021e72a00SBarry Smith 
65121e72a00SBarry Smith    Input Parameters:
65221e72a00SBarry Smith +  an  - number of values in the first array
65321e72a00SBarry Smith .  aI  - first sorted array of integers
65421e72a00SBarry Smith .  bn  - number of values in the second array
65521e72a00SBarry Smith -  bI  - second array of integers
65621e72a00SBarry Smith 
65721e72a00SBarry Smith    Output Parameters:
65821e72a00SBarry Smith +  n   - number of values in the merged array
65921e72a00SBarry Smith -  I   - merged sorted array, this is allocated if an array is not provided
66021e72a00SBarry Smith 
66121e72a00SBarry Smith    Level: intermediate
66221e72a00SBarry Smith 
66321e72a00SBarry Smith    Concepts: merging^arrays
66421e72a00SBarry Smith 
66521e72a00SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
66621e72a00SBarry Smith @*/
66721e72a00SBarry Smith PetscErrorCode  PetscMergeIntArray(PetscInt an,const PetscInt *aI, PetscInt bn, const PetscInt *bI,  PetscInt *n, PetscInt **L)
66821e72a00SBarry Smith {
66921e72a00SBarry Smith   PetscErrorCode ierr;
67021e72a00SBarry Smith   PetscInt       *L_ = *L, ak, bk, k;
67121e72a00SBarry Smith 
67221e72a00SBarry Smith   if (!L_) {
67321e72a00SBarry Smith     ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr);
67421e72a00SBarry Smith     L_   = *L;
67521e72a00SBarry Smith   }
67621e72a00SBarry Smith   k = ak = bk = 0;
67721e72a00SBarry Smith   while (ak < an && bk < bn) {
67821e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
67921e72a00SBarry Smith       L_[k] = aI[ak];
68021e72a00SBarry Smith       ++ak;
68121e72a00SBarry Smith       ++bk;
68221e72a00SBarry Smith       ++k;
68321e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
68421e72a00SBarry Smith       L_[k] = aI[ak];
68521e72a00SBarry Smith       ++ak;
68621e72a00SBarry Smith       ++k;
68721e72a00SBarry Smith     } else {
68821e72a00SBarry Smith       L_[k] = bI[bk];
68921e72a00SBarry Smith       ++bk;
69021e72a00SBarry Smith       ++k;
69121e72a00SBarry Smith     }
69221e72a00SBarry Smith   }
69321e72a00SBarry Smith   if (ak < an) {
69421e72a00SBarry Smith     ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr);
69521e72a00SBarry Smith     k   += (an-ak);
69621e72a00SBarry Smith   }
69721e72a00SBarry Smith   if (bk < bn) {
69821e72a00SBarry Smith     ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr);
69921e72a00SBarry Smith     k   += (bn-bk);
70021e72a00SBarry Smith   }
70121e72a00SBarry Smith   *n = k;
70221e72a00SBarry Smith   PetscFunctionReturn(0);
70321e72a00SBarry Smith }
704b4301105SBarry Smith 
705b4301105SBarry Smith #undef __FUNCT__
706c1f0200aSDmitry Karpeev #define __FUNCT__ "PetscMergeIntArrayPair"
707c1f0200aSDmitry Karpeev /*@
70821e72a00SBarry Smith    PetscMergeIntArrayPair -     Merges two SORTED integer arrays that share NO common values along with an additional array of integers.
709c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
710c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
711c1f0200aSDmitry Karpeev 
712c1f0200aSDmitry Karpeev    Not Collective
713c1f0200aSDmitry Karpeev 
714c1f0200aSDmitry Karpeev    Input Parameters:
715c1f0200aSDmitry Karpeev +  an  - number of values in the first array
716c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
717c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
718c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
719c1f0200aSDmitry Karpeev .  bI  - second array of integers
720c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
721c1f0200aSDmitry Karpeev 
722c1f0200aSDmitry Karpeev    Output Parameters:
723c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
724c1f0200aSDmitry Karpeev .  I   - merged sorted array
725c1f0200aSDmitry Karpeev -  J   - merged additional array
726c1f0200aSDmitry Karpeev 
727c1f0200aSDmitry Karpeev    Level: intermediate
728c1f0200aSDmitry Karpeev 
729c1f0200aSDmitry Karpeev    Concepts: merging^arrays
730c1f0200aSDmitry Karpeev 
731c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
732c1f0200aSDmitry Karpeev @*/
733d7aa01a8SHong 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)
734c1f0200aSDmitry Karpeev {
735c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
736d7aa01a8SHong Zhang   PetscInt       n_, *L_ = *L, *J_= *J, ak, bk, k;
737c1f0200aSDmitry Karpeev 
738c1f0200aSDmitry Karpeev   n_ = an + bn;
739c1f0200aSDmitry Karpeev   *n = n_;
740d7aa01a8SHong Zhang   if (!L_) {
741785e854fSJed Brown     ierr = PetscMalloc1(n_, L);CHKERRQ(ierr);
742d7aa01a8SHong Zhang     L_   = *L;
743c1f0200aSDmitry Karpeev   }
744c1f0200aSDmitry Karpeev   if (!J_) {
745785e854fSJed Brown     ierr = PetscMalloc1(n_, &J_);CHKERRQ(ierr);
746c1f0200aSDmitry Karpeev     J_   = *J;
747c1f0200aSDmitry Karpeev   }
748c1f0200aSDmitry Karpeev   k = ak = bk = 0;
749c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
750c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
751d7aa01a8SHong Zhang       L_[k] = aI[ak];
752c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
753c1f0200aSDmitry Karpeev       ++ak;
754c1f0200aSDmitry Karpeev       ++k;
755a297a907SKarl Rupp     } else {
756d7aa01a8SHong Zhang       L_[k] = bI[bk];
757c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
758c1f0200aSDmitry Karpeev       ++bk;
759c1f0200aSDmitry Karpeev       ++k;
760c1f0200aSDmitry Karpeev     }
761c1f0200aSDmitry Karpeev   }
762c1f0200aSDmitry Karpeev   if (ak < an) {
763d7aa01a8SHong Zhang     ierr = PetscMemcpy(L_+k,aI+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr);
764c1f0200aSDmitry Karpeev     ierr = PetscMemcpy(J_+k,aJ+ak,(an-ak)*sizeof(PetscInt));CHKERRQ(ierr);
765c1f0200aSDmitry Karpeev     k   += (an-ak);
766c1f0200aSDmitry Karpeev   }
767c1f0200aSDmitry Karpeev   if (bk < bn) {
768d7aa01a8SHong Zhang     ierr = PetscMemcpy(L_+k,bI+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr);
769c1f0200aSDmitry Karpeev     ierr = PetscMemcpy(J_+k,bJ+bk,(bn-bk)*sizeof(PetscInt));CHKERRQ(ierr);
770c1f0200aSDmitry Karpeev   }
771c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
772c1f0200aSDmitry Karpeev }
773c1f0200aSDmitry Karpeev 
774e5c89e4eSSatish Balay 
77542eaa7f3SBarry Smith #undef __FUNCT__
77642eaa7f3SBarry Smith #define __FUNCT__ "PetscProcessTree"
77742eaa7f3SBarry Smith /*@
77842eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
77942eaa7f3SBarry Smith 
78042eaa7f3SBarry Smith    Not Collective
78142eaa7f3SBarry Smith 
78242eaa7f3SBarry Smith    Input Parameters:
78342eaa7f3SBarry Smith +  n  - number of values
78442eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
78542eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
78642eaa7f3SBarry Smith 
78742eaa7f3SBarry Smith    Output Parameters:
78842eaa7f3SBarry Smith +  Nlevels - the number of levels
78942eaa7f3SBarry Smith .  Level - for each node tells its level
79042eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
79142eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
79242eaa7f3SBarry Smith -  Column - for each id tells its column index
79342eaa7f3SBarry Smith 
79442eaa7f3SBarry Smith    Level: intermediate
79542eaa7f3SBarry Smith 
79642eaa7f3SBarry Smith 
79742eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation()
79842eaa7f3SBarry Smith @*/
7997087cfbeSBarry Smith PetscErrorCode  PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column)
80042eaa7f3SBarry Smith {
80142eaa7f3SBarry Smith   PetscInt       i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column;
80242eaa7f3SBarry Smith   PetscErrorCode ierr;
803ace3abfcSBarry Smith   PetscBool      done = PETSC_FALSE;
80442eaa7f3SBarry Smith 
80542eaa7f3SBarry Smith   PetscFunctionBegin;
80642eaa7f3SBarry Smith   if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set");
80742eaa7f3SBarry Smith   for (i=0; i<n; i++) {
80842eaa7f3SBarry Smith     if (mask[i]) continue;
80942eaa7f3SBarry Smith     if (parentid[i]  == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent");
81042eaa7f3SBarry Smith     if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked");
81142eaa7f3SBarry Smith   }
81242eaa7f3SBarry Smith 
81342eaa7f3SBarry Smith   for (i=0; i<n; i++) {
81442eaa7f3SBarry Smith     if (!mask[i]) nmask++;
81542eaa7f3SBarry Smith   }
81642eaa7f3SBarry Smith 
81742eaa7f3SBarry Smith   /* determine the level in the tree of each node */
8181795a4d1SJed Brown   ierr = PetscCalloc1(n,&level);CHKERRQ(ierr);
819a297a907SKarl Rupp 
82042eaa7f3SBarry Smith   level[0] = 1;
82142eaa7f3SBarry Smith   while (!done) {
82242eaa7f3SBarry Smith     done = PETSC_TRUE;
82342eaa7f3SBarry Smith     for (i=0; i<n; i++) {
82442eaa7f3SBarry Smith       if (mask[i]) continue;
82542eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
82642eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
82742eaa7f3SBarry Smith     }
82842eaa7f3SBarry Smith   }
82942eaa7f3SBarry Smith   for (i=0; i<n; i++) {
83042eaa7f3SBarry Smith     level[i]--;
83142eaa7f3SBarry Smith     nlevels = PetscMax(nlevels,level[i]);
83242eaa7f3SBarry Smith   }
83342eaa7f3SBarry Smith 
83442eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
8351795a4d1SJed Brown   ierr = PetscCalloc1(nlevels,&levelcnt);CHKERRQ(ierr);
83642eaa7f3SBarry Smith   for (i=0; i<n; i++) {
83742eaa7f3SBarry Smith     if (mask[i]) continue;
83842eaa7f3SBarry Smith     levelcnt[level[i]-1]++;
83942eaa7f3SBarry Smith   }
840a297a907SKarl Rupp   for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]);
84142eaa7f3SBarry Smith 
84242eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
843dcca6d9dSJed Brown   ierr = PetscMalloc2(levelmax,&workid,levelmax,&workparentid);CHKERRQ(ierr);
844785e854fSJed Brown   ierr = PetscMalloc1(nmask,&idbylevel);CHKERRQ(ierr);
84542eaa7f3SBarry Smith   for (j=1; j<=nlevels;j++) {
84642eaa7f3SBarry Smith     cnt = 0;
84742eaa7f3SBarry Smith     for (i=0; i<n; i++) {
84842eaa7f3SBarry Smith       if (mask[i]) continue;
84942eaa7f3SBarry Smith       if (level[i] != j) continue;
85042eaa7f3SBarry Smith       workid[cnt]         = i;
85142eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
85242eaa7f3SBarry Smith     }
85342eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
85442eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
85542eaa7f3SBarry Smith     ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr);
85642eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
85742eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
85842eaa7f3SBarry Smith     ierr  = PetscMemcpy(idbylevel+tcnt,workid,cnt*sizeof(PetscInt));CHKERRQ(ierr);
85942eaa7f3SBarry Smith     tcnt += cnt;
86042eaa7f3SBarry Smith   }
86142eaa7f3SBarry Smith   if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes");
86242eaa7f3SBarry Smith   ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr);
86342eaa7f3SBarry Smith 
86442eaa7f3SBarry Smith   /* for each node list its column */
865785e854fSJed Brown   ierr = PetscMalloc1(n,&column);CHKERRQ(ierr);
86642eaa7f3SBarry Smith   cnt = 0;
86742eaa7f3SBarry Smith   for (j=0; j<nlevels; j++) {
86842eaa7f3SBarry Smith     for (i=0; i<levelcnt[j]; i++) {
86942eaa7f3SBarry Smith       column[idbylevel[cnt++]] = i;
87042eaa7f3SBarry Smith     }
87142eaa7f3SBarry Smith   }
87242eaa7f3SBarry Smith 
87342eaa7f3SBarry Smith   *Nlevels   = nlevels;
87442eaa7f3SBarry Smith   *Level     = level;
87542eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
87642eaa7f3SBarry Smith   *Idbylevel = idbylevel;
87742eaa7f3SBarry Smith   *Column    = column;
87842eaa7f3SBarry Smith   PetscFunctionReturn(0);
87942eaa7f3SBarry Smith }
880