xref: /petsc/src/sys/utils/sorti.c (revision 59e16d97b8ec1d363b4b1a34d3190725b7824d5a)
17d0a6c19SBarry Smith 
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file contains routines for sorting integers. Values are sorted in place.
4*59e16d97SJunchao Zhang    One can use src/sys/examples/tests/ex52.c for benchmarking.
5e5c89e4eSSatish Balay  */
6af0996ceSBarry Smith #include <petsc/private/petscimpl.h>                /*I  "petscsys.h"  I*/
7e5c89e4eSSatish Balay 
8a8582168SJed Brown #define MEDIAN3(v,a,b,c)                                                        \
9a8582168SJed Brown   (v[a]<v[b]                                                                    \
10a8582168SJed Brown    ? (v[b]<v[c]                                                                 \
11*59e16d97SJunchao Zhang       ? (b)                                                                     \
12*59e16d97SJunchao Zhang       : (v[a]<v[c] ? (c) : (a)))                                                \
13a8582168SJed Brown    : (v[c]<v[b]                                                                 \
14*59e16d97SJunchao Zhang       ? (b)                                                                     \
15*59e16d97SJunchao Zhang       : (v[a]<v[c] ? (a) : (c))))
16a8582168SJed Brown 
17a8582168SJed Brown #define MEDIAN(v,right) MEDIAN3(v,right/4,right/2,right/4*3)
18ef8e3583SJed Brown 
19*59e16d97SJunchao Zhang /* Swap one, two or three pairs. Each pair can have its own type */
20*59e16d97SJunchao Zhang #define SWAP1(a,b,t1)               do {t1=a;a=b;b=t1;} while(0)
21*59e16d97SJunchao Zhang #define SWAP2(a,b,c,d,t1,t2)        do {t1=a;a=b;b=t1; t2=c;c=d;d=t2;} while(0)
22*59e16d97SJunchao Zhang #define SWAP3(a,b,c,d,e,f,t1,t2,t3) do {t1=a;a=b;b=t1; t2=c;c=d;d=t2; t3=e;e=f;f=t3;} while(0)
23*59e16d97SJunchao Zhang 
24*59e16d97SJunchao Zhang /* Swap a & b, *c & *d. c, d, t2 are pointers to a type of size <siz> */
25*59e16d97SJunchao Zhang #define SWAP2Data(a,b,c,d,t1,t2,siz)                                             \
26*59e16d97SJunchao Zhang   do {                                                                           \
27*59e16d97SJunchao Zhang     PetscErrorCode ierr;                                                         \
28*59e16d97SJunchao Zhang     t1=a; a=b; b=t1;                                                             \
29*59e16d97SJunchao Zhang     ierr = PetscMemcpy(t2,c,siz);CHKERRQ(ierr);                                  \
30*59e16d97SJunchao Zhang     ierr = PetscMemcpy(c,d,siz);CHKERRQ(ierr);                                   \
31*59e16d97SJunchao Zhang     ierr = PetscMemcpy(d,t2,siz);CHKERRQ(ierr);                                  \
32*59e16d97SJunchao Zhang   } while(0)
33e5c89e4eSSatish Balay 
34e5c89e4eSSatish Balay /*
35*59e16d97SJunchao Zhang    Partition X[lo,hi] into two parts: X[lo,l) <= pivot; X[r,hi] > pivot
36e5c89e4eSSatish Balay 
37*59e16d97SJunchao Zhang    Input Parameters:
38*59e16d97SJunchao Zhang     + X         - array to partition
39*59e16d97SJunchao Zhang     . pivot     - a pivot of X[]
40*59e16d97SJunchao Zhang     . t1        - temp variable for X
41*59e16d97SJunchao Zhang     - lo,hi     - lower and upper bound of the array
42*59e16d97SJunchao Zhang 
43*59e16d97SJunchao Zhang    Output Parameters:
44*59e16d97SJunchao Zhang     + l,r       - of type PetscInt
45*59e16d97SJunchao Zhang 
46*59e16d97SJunchao Zhang    Notes:
47*59e16d97SJunchao Zhang     The TwoWayPartition2/3 variants also partition other arrays along with X.
48*59e16d97SJunchao Zhang     These arrays can have different types, so they provide their own temp t2,t3
49*59e16d97SJunchao Zhang  */
50*59e16d97SJunchao Zhang #define TwoWayPartition1(X,pivot,t1,lo,hi,l,r)                                   \
51*59e16d97SJunchao Zhang   do {                                                                           \
52*59e16d97SJunchao Zhang     l = lo;                                                                      \
53*59e16d97SJunchao Zhang     r = hi;                                                                      \
54*59e16d97SJunchao Zhang     while(1) {                                                                   \
55*59e16d97SJunchao Zhang       while (X[l] < pivot) l++;                                                  \
56*59e16d97SJunchao Zhang       while (X[r] > pivot) r--;                                                  \
57*59e16d97SJunchao Zhang       if (l >= r) {r++; break;}                                                  \
58*59e16d97SJunchao Zhang       SWAP1(X[l],X[r],t1);                                                       \
59*59e16d97SJunchao Zhang       l++;                                                                       \
60*59e16d97SJunchao Zhang       r--;                                                                       \
61*59e16d97SJunchao Zhang     }                                                                            \
62*59e16d97SJunchao Zhang   } while(0)
63*59e16d97SJunchao Zhang 
64*59e16d97SJunchao Zhang #define TwoWayPartition2(X,Y,pivot,t1,t2,lo,hi,l,r)                              \
65*59e16d97SJunchao Zhang   do {                                                                           \
66*59e16d97SJunchao Zhang     l = lo;                                                                      \
67*59e16d97SJunchao Zhang     r = hi;                                                                      \
68*59e16d97SJunchao Zhang     while(1) {                                                                   \
69*59e16d97SJunchao Zhang       while (X[l] < pivot) l++;                                                  \
70*59e16d97SJunchao Zhang       while (X[r] > pivot) r--;                                                  \
71*59e16d97SJunchao Zhang       if (l >= r) {r++; break;}                                                  \
72*59e16d97SJunchao Zhang       SWAP2(X[l],X[r],Y[l],Y[r],t1,t2);                                          \
73*59e16d97SJunchao Zhang       l++;                                                                       \
74*59e16d97SJunchao Zhang       r--;                                                                       \
75*59e16d97SJunchao Zhang     }                                                                            \
76*59e16d97SJunchao Zhang   } while(0)
77*59e16d97SJunchao Zhang 
78*59e16d97SJunchao Zhang #define TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,lo,hi,l,r)                         \
79*59e16d97SJunchao Zhang   do {                                                                           \
80*59e16d97SJunchao Zhang     l = lo;                                                                      \
81*59e16d97SJunchao Zhang     r = hi;                                                                      \
82*59e16d97SJunchao Zhang     while(1) {                                                                   \
83*59e16d97SJunchao Zhang       while (X[l] < pivot) l++;                                                  \
84*59e16d97SJunchao Zhang       while (X[r] > pivot) r--;                                                  \
85*59e16d97SJunchao Zhang       if (l >= r) {r++; break;}                                                  \
86*59e16d97SJunchao Zhang       SWAP3(X[l],X[r],Y[l],Y[r],Z[l],Z[r],t1,t2,t3);                             \
87*59e16d97SJunchao Zhang       l++;                                                                       \
88*59e16d97SJunchao Zhang       r--;                                                                       \
89*59e16d97SJunchao Zhang     }                                                                            \
90*59e16d97SJunchao Zhang   } while(0)
91*59e16d97SJunchao Zhang 
92*59e16d97SJunchao Zhang /* Templates for similar functions used below */
93*59e16d97SJunchao Zhang #define QuickSort1(FuncName,X,n,pivot,t1,ierr)                                   \
94*59e16d97SJunchao Zhang   do {                                                                           \
95*59e16d97SJunchao Zhang     PetscInt i,j,p,l,r,hi=n-1;                                                   \
96*59e16d97SJunchao Zhang     if (n < 8) {                                                                 \
97*59e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
98*59e16d97SJunchao Zhang         pivot = X[i];                                                            \
99*59e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
100*59e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
101*59e16d97SJunchao Zhang             SWAP1(X[i],X[j],t1);                                                 \
102*59e16d97SJunchao Zhang             pivot = X[i];                                                        \
103*59e16d97SJunchao Zhang           }                                                                      \
104*59e16d97SJunchao Zhang         }                                                                        \
105*59e16d97SJunchao Zhang       }                                                                          \
106*59e16d97SJunchao Zhang     } else {                                                                     \
107*59e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
108*59e16d97SJunchao Zhang       pivot = X[p];                                                              \
109*59e16d97SJunchao Zhang       TwoWayPartition1(X,pivot,t1,0,hi,l,r);                                     \
110*59e16d97SJunchao Zhang       ierr  = FuncName(l,X);CHKERRQ(ierr);                                       \
111*59e16d97SJunchao Zhang       ierr  = FuncName(hi-r+1,X+r);CHKERRQ(ierr);                                \
112*59e16d97SJunchao Zhang     }                                                                            \
113*59e16d97SJunchao Zhang   } while(0)
114*59e16d97SJunchao Zhang 
115*59e16d97SJunchao Zhang #define QuickSort2(FuncName,X,Y,n,pivot,t1,t2,ierr)                              \
116*59e16d97SJunchao Zhang   do {                                                                           \
117*59e16d97SJunchao Zhang     PetscInt i,j,p,l,r,hi=n-1;                                                   \
118*59e16d97SJunchao Zhang     if (n < 8) {                                                                 \
119*59e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
120*59e16d97SJunchao Zhang         pivot = X[i];                                                            \
121*59e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
122*59e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
123*59e16d97SJunchao Zhang             SWAP2(X[i],X[j],Y[i],Y[j],t1,t2);                                    \
124*59e16d97SJunchao Zhang             pivot = X[i];                                                        \
125*59e16d97SJunchao Zhang           }                                                                      \
126*59e16d97SJunchao Zhang         }                                                                        \
127*59e16d97SJunchao Zhang       }                                                                          \
128*59e16d97SJunchao Zhang     } else {                                                                     \
129*59e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
130*59e16d97SJunchao Zhang       pivot = X[p];                                                              \
131*59e16d97SJunchao Zhang       TwoWayPartition2(X,Y,pivot,t1,t2,0,hi,l,r);                                \
132*59e16d97SJunchao Zhang       ierr  = FuncName(l,X,Y);CHKERRQ(ierr);                                     \
133*59e16d97SJunchao Zhang       ierr  = FuncName(hi-r+1,X+r,Y+r);CHKERRQ(ierr);                            \
134*59e16d97SJunchao Zhang     }                                                                            \
135*59e16d97SJunchao Zhang   } while(0)
136*59e16d97SJunchao Zhang 
137*59e16d97SJunchao Zhang #define QuickSort3(FuncName,X,Y,Z,n,pivot,t1,t2,t3,ierr)                         \
138*59e16d97SJunchao Zhang   do {                                                                           \
139*59e16d97SJunchao Zhang     PetscInt i,j,p,l,r,hi=n-1;                                                   \
140*59e16d97SJunchao Zhang     if (n < 8) {                                                                 \
141*59e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
142*59e16d97SJunchao Zhang         pivot = X[i];                                                            \
143*59e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
144*59e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
145*59e16d97SJunchao Zhang             SWAP3(X[i],X[j],Y[i],Y[j],Z[i],Z[j],t1,t2,t3);                       \
146*59e16d97SJunchao Zhang             pivot = X[i];                                                        \
147*59e16d97SJunchao Zhang           }                                                                      \
148*59e16d97SJunchao Zhang         }                                                                        \
149*59e16d97SJunchao Zhang       }                                                                          \
150*59e16d97SJunchao Zhang     } else {                                                                     \
151*59e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
152*59e16d97SJunchao Zhang       pivot = X[p];                                                              \
153*59e16d97SJunchao Zhang       TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,0,hi,l,r);                           \
154*59e16d97SJunchao Zhang       ierr  = FuncName(l,X,Y,Z);CHKERRQ(ierr);                                   \
155*59e16d97SJunchao Zhang       ierr  = FuncName(hi-r+1,X+r,Y+r,Z+r);CHKERRQ(ierr);                        \
156*59e16d97SJunchao Zhang     }                                                                            \
157*59e16d97SJunchao Zhang   } while(0)
158e5c89e4eSSatish Balay 
159e5c89e4eSSatish Balay /*@
160e5c89e4eSSatish Balay    PetscSortInt - Sorts an array of integers in place in increasing order.
161e5c89e4eSSatish Balay 
162e5c89e4eSSatish Balay    Not Collective
163e5c89e4eSSatish Balay 
164e5c89e4eSSatish Balay    Input Parameters:
165e5c89e4eSSatish Balay +  n  - number of values
166*59e16d97SJunchao Zhang -  X  - array of integers
167e5c89e4eSSatish Balay 
168e5c89e4eSSatish Balay    Level: intermediate
169e5c89e4eSSatish Balay 
170e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntWithPermutation()
171e5c89e4eSSatish Balay @*/
172*59e16d97SJunchao Zhang PetscErrorCode  PetscSortInt(PetscInt n,PetscInt *X)
173e5c89e4eSSatish Balay {
174*59e16d97SJunchao Zhang   PetscErrorCode ierr;
175*59e16d97SJunchao Zhang   PetscInt       pivot,t1;
176e5c89e4eSSatish Balay 
177e5c89e4eSSatish Balay   PetscFunctionBegin;
178*59e16d97SJunchao Zhang   QuickSort1(PetscSortInt,X,n,pivot,t1,ierr);
179e5c89e4eSSatish Balay   PetscFunctionReturn(0);
180e5c89e4eSSatish Balay }
181e5c89e4eSSatish Balay 
1821db36a52SBarry Smith /*@
18322ab5688SLisandro Dalcin    PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted input array
18422ab5688SLisandro Dalcin 
18522ab5688SLisandro Dalcin    Not Collective
18622ab5688SLisandro Dalcin 
18722ab5688SLisandro Dalcin    Input Parameters:
18822ab5688SLisandro Dalcin +  n  - number of values
189*59e16d97SJunchao Zhang -  X  - sorted array of integers
19022ab5688SLisandro Dalcin 
19122ab5688SLisandro Dalcin    Output Parameter:
19222ab5688SLisandro Dalcin .  n - number of non-redundant values
19322ab5688SLisandro Dalcin 
19422ab5688SLisandro Dalcin    Level: intermediate
19522ab5688SLisandro Dalcin 
19622ab5688SLisandro Dalcin .seealso: PetscSortInt()
19722ab5688SLisandro Dalcin @*/
198*59e16d97SJunchao Zhang PetscErrorCode  PetscSortedRemoveDupsInt(PetscInt *n,PetscInt X[])
19922ab5688SLisandro Dalcin {
20022ab5688SLisandro Dalcin   PetscInt i,s = 0,N = *n, b = 0;
20122ab5688SLisandro Dalcin 
20222ab5688SLisandro Dalcin   PetscFunctionBegin;
20322ab5688SLisandro Dalcin   for (i=0; i<N-1; i++) {
204*59e16d97SJunchao Zhang     if (PetscUnlikely(X[b+s+1] < X[b])) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Input array is not sorted");
205*59e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
206*59e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
20722ab5688SLisandro Dalcin     } else s++;
20822ab5688SLisandro Dalcin   }
20922ab5688SLisandro Dalcin   *n = N - s;
21022ab5688SLisandro Dalcin   PetscFunctionReturn(0);
21122ab5688SLisandro Dalcin }
21222ab5688SLisandro Dalcin 
21322ab5688SLisandro Dalcin /*@
2141db36a52SBarry Smith    PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries
2151db36a52SBarry Smith 
2161db36a52SBarry Smith    Not Collective
2171db36a52SBarry Smith 
2181db36a52SBarry Smith    Input Parameters:
2191db36a52SBarry Smith +  n  - number of values
220*59e16d97SJunchao Zhang -  X  - array of integers
2211db36a52SBarry Smith 
2221db36a52SBarry Smith    Output Parameter:
2231db36a52SBarry Smith .  n - number of non-redundant values
2241db36a52SBarry Smith 
2251db36a52SBarry Smith    Level: intermediate
2261db36a52SBarry Smith 
22722ab5688SLisandro Dalcin .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortedRemoveDupsInt()
2281db36a52SBarry Smith @*/
229*59e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsInt(PetscInt *n,PetscInt X[])
2301db36a52SBarry Smith {
2311db36a52SBarry Smith   PetscErrorCode ierr;
2321db36a52SBarry Smith 
2331db36a52SBarry Smith   PetscFunctionBegin;
234*59e16d97SJunchao Zhang   ierr = PetscSortInt(*n,X);CHKERRQ(ierr);
235*59e16d97SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(n,X);CHKERRQ(ierr);
2361db36a52SBarry Smith   PetscFunctionReturn(0);
2371db36a52SBarry Smith }
2381db36a52SBarry Smith 
23960e03357SMatthew G Knepley /*@
240d9f1162dSJed Brown   PetscFindInt - Finds integer in a sorted array of integers
24160e03357SMatthew G Knepley 
24260e03357SMatthew G Knepley    Not Collective
24360e03357SMatthew G Knepley 
24460e03357SMatthew G Knepley    Input Parameters:
24560e03357SMatthew G Knepley +  key - the integer to locate
246d9f1162dSJed Brown .  n   - number of values in the array
247*59e16d97SJunchao Zhang -  X  - array of integers
24860e03357SMatthew G Knepley 
24960e03357SMatthew G Knepley    Output Parameter:
250d9f1162dSJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
25160e03357SMatthew G Knepley 
25260e03357SMatthew G Knepley    Level: intermediate
25360e03357SMatthew G Knepley 
25460e03357SMatthew G Knepley .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
25560e03357SMatthew G Knepley @*/
256*59e16d97SJunchao Zhang PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
25760e03357SMatthew G Knepley {
258d9f1162dSJed Brown   PetscInt lo = 0,hi = n;
25960e03357SMatthew G Knepley 
26060e03357SMatthew G Knepley   PetscFunctionBegin;
26160e03357SMatthew G Knepley   PetscValidPointer(loc,4);
26298781d82SMatthew G Knepley   if (!n) {*loc = -1; PetscFunctionReturn(0);}
263*59e16d97SJunchao Zhang   PetscValidPointer(X,3);
264d9f1162dSJed Brown   while (hi - lo > 1) {
265d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo)/2;
266*59e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
267d9f1162dSJed Brown     else               lo = mid;
26860e03357SMatthew G Knepley   }
269*59e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
27060e03357SMatthew G Knepley   PetscFunctionReturn(0);
27160e03357SMatthew G Knepley }
27260e03357SMatthew G Knepley 
273d2aeb606SJed Brown /*@
274d2aeb606SJed Brown   PetscFindMPIInt - Finds MPI integer in a sorted array of integers
275d2aeb606SJed Brown 
276d2aeb606SJed Brown    Not Collective
277d2aeb606SJed Brown 
278d2aeb606SJed Brown    Input Parameters:
279d2aeb606SJed Brown +  key - the integer to locate
280d2aeb606SJed Brown .  n   - number of values in the array
281*59e16d97SJunchao Zhang -  X   - array of integers
282d2aeb606SJed Brown 
283d2aeb606SJed Brown    Output Parameter:
284d2aeb606SJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
285d2aeb606SJed Brown 
286d2aeb606SJed Brown    Level: intermediate
287d2aeb606SJed Brown 
288d2aeb606SJed Brown .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
289d2aeb606SJed Brown @*/
290*59e16d97SJunchao Zhang PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
291d2aeb606SJed Brown {
292d2aeb606SJed Brown   PetscInt lo = 0,hi = n;
293d2aeb606SJed Brown 
294d2aeb606SJed Brown   PetscFunctionBegin;
295d2aeb606SJed Brown   PetscValidPointer(loc,4);
296d2aeb606SJed Brown   if (!n) {*loc = -1; PetscFunctionReturn(0);}
297*59e16d97SJunchao Zhang   PetscValidPointer(X,3);
298d2aeb606SJed Brown   while (hi - lo > 1) {
299d2aeb606SJed Brown     PetscInt mid = lo + (hi - lo)/2;
300*59e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
301d2aeb606SJed Brown     else               lo = mid;
302d2aeb606SJed Brown   }
303*59e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
304e5c89e4eSSatish Balay   PetscFunctionReturn(0);
305e5c89e4eSSatish Balay }
306e5c89e4eSSatish Balay 
307e5c89e4eSSatish Balay /*@
308e5c89e4eSSatish Balay    PetscSortIntWithArray - Sorts an array of integers in place in increasing order;
309e5c89e4eSSatish Balay        changes a second array to match the sorted first array.
310e5c89e4eSSatish Balay 
311e5c89e4eSSatish Balay    Not Collective
312e5c89e4eSSatish Balay 
313e5c89e4eSSatish Balay    Input Parameters:
314e5c89e4eSSatish Balay +  n  - number of values
315*59e16d97SJunchao Zhang .  X  - array of integers
316*59e16d97SJunchao Zhang -  Y  - second array of integers
317e5c89e4eSSatish Balay 
318e5c89e4eSSatish Balay    Level: intermediate
319e5c89e4eSSatish Balay 
320e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
321e5c89e4eSSatish Balay @*/
322*59e16d97SJunchao Zhang PetscErrorCode  PetscSortIntWithArray(PetscInt n,PetscInt *X,PetscInt *Y)
323e5c89e4eSSatish Balay {
324e5c89e4eSSatish Balay   PetscErrorCode ierr;
325*59e16d97SJunchao Zhang   PetscInt       pivot,t1,t2;
326e5c89e4eSSatish Balay 
327e5c89e4eSSatish Balay   PetscFunctionBegin;
328*59e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithArray,X,Y,n,pivot,t1,t2,ierr);
329c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
330c1f0200aSDmitry Karpeev }
331c1f0200aSDmitry Karpeev 
332c1f0200aSDmitry Karpeev /*@
333c1f0200aSDmitry Karpeev    PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order;
334c1f0200aSDmitry Karpeev        changes a pair of integer arrays to match the sorted first array.
335c1f0200aSDmitry Karpeev 
336c1f0200aSDmitry Karpeev    Not Collective
337c1f0200aSDmitry Karpeev 
338c1f0200aSDmitry Karpeev    Input Parameters:
339c1f0200aSDmitry Karpeev +  n  - number of values
340*59e16d97SJunchao Zhang .  X  - array of integers
341*59e16d97SJunchao Zhang .  Y  - second array of integers (first array of the pair)
342*59e16d97SJunchao Zhang -  Z  - third array of integers  (second array of the pair)
343c1f0200aSDmitry Karpeev 
344c1f0200aSDmitry Karpeev    Level: intermediate
345c1f0200aSDmitry Karpeev 
346c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray()
347c1f0200aSDmitry Karpeev @*/
348*59e16d97SJunchao Zhang PetscErrorCode  PetscSortIntWithArrayPair(PetscInt n,PetscInt *X,PetscInt *Y,PetscInt *Z)
349c1f0200aSDmitry Karpeev {
350c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
351*59e16d97SJunchao Zhang   PetscInt       pivot,t1,t2,t3;
352c1f0200aSDmitry Karpeev 
353c1f0200aSDmitry Karpeev   PetscFunctionBegin;
354*59e16d97SJunchao Zhang   QuickSort3(PetscSortIntWithArrayPair,X,Y,Z,n,pivot,t1,t2,t3,ierr);
355c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
356c1f0200aSDmitry Karpeev }
357c1f0200aSDmitry Karpeev 
35817d7d925SStefano Zampini /*@
35917d7d925SStefano Zampini    PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order.
36017d7d925SStefano Zampini 
36117d7d925SStefano Zampini    Not Collective
36217d7d925SStefano Zampini 
36317d7d925SStefano Zampini    Input Parameters:
36417d7d925SStefano Zampini +  n  - number of values
365*59e16d97SJunchao Zhang -  X  - array of integers
36617d7d925SStefano Zampini 
36717d7d925SStefano Zampini    Level: intermediate
36817d7d925SStefano Zampini 
36917d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation()
37017d7d925SStefano Zampini @*/
371*59e16d97SJunchao Zhang PetscErrorCode  PetscSortMPIInt(PetscInt n,PetscMPIInt *X)
37217d7d925SStefano Zampini {
373*59e16d97SJunchao Zhang   PetscErrorCode ierr;
374*59e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
37517d7d925SStefano Zampini 
37617d7d925SStefano Zampini   PetscFunctionBegin;
377*59e16d97SJunchao Zhang   QuickSort1(PetscSortMPIInt,X,n,pivot,t1,ierr);
37817d7d925SStefano Zampini   PetscFunctionReturn(0);
37917d7d925SStefano Zampini }
38017d7d925SStefano Zampini 
38117d7d925SStefano Zampini /*@
38217d7d925SStefano Zampini    PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries
38317d7d925SStefano Zampini 
38417d7d925SStefano Zampini    Not Collective
38517d7d925SStefano Zampini 
38617d7d925SStefano Zampini    Input Parameters:
38717d7d925SStefano Zampini +  n  - number of values
388*59e16d97SJunchao Zhang -  X  - array of integers
38917d7d925SStefano Zampini 
39017d7d925SStefano Zampini    Output Parameter:
39117d7d925SStefano Zampini .  n - number of non-redundant values
39217d7d925SStefano Zampini 
39317d7d925SStefano Zampini    Level: intermediate
39417d7d925SStefano Zampini 
39517d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
39617d7d925SStefano Zampini @*/
397*59e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt X[])
39817d7d925SStefano Zampini {
39917d7d925SStefano Zampini   PetscErrorCode ierr;
40017d7d925SStefano Zampini   PetscInt       i,s = 0,N = *n, b = 0;
40117d7d925SStefano Zampini 
40217d7d925SStefano Zampini   PetscFunctionBegin;
403*59e16d97SJunchao Zhang   ierr = PetscSortMPIInt(N,X);CHKERRQ(ierr);
40417d7d925SStefano Zampini   for (i=0; i<N-1; i++) {
405*59e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
406*59e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
407a297a907SKarl Rupp     } else s++;
40817d7d925SStefano Zampini   }
40917d7d925SStefano Zampini   *n = N - s;
41017d7d925SStefano Zampini   PetscFunctionReturn(0);
41117d7d925SStefano Zampini }
412c1f0200aSDmitry Karpeev 
4134d615ea0SBarry Smith /*@
4144d615ea0SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order;
4154d615ea0SBarry Smith        changes a second array to match the sorted first array.
4164d615ea0SBarry Smith 
4174d615ea0SBarry Smith    Not Collective
4184d615ea0SBarry Smith 
4194d615ea0SBarry Smith    Input Parameters:
4204d615ea0SBarry Smith +  n  - number of values
421*59e16d97SJunchao Zhang .  X  - array of integers
422*59e16d97SJunchao Zhang -  Y  - second array of integers
4234d615ea0SBarry Smith 
4244d615ea0SBarry Smith    Level: intermediate
4254d615ea0SBarry Smith 
4264d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
4274d615ea0SBarry Smith @*/
428*59e16d97SJunchao Zhang PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt *X,PetscMPIInt *Y)
4294d615ea0SBarry Smith {
4304d615ea0SBarry Smith   PetscErrorCode ierr;
431*59e16d97SJunchao Zhang   PetscMPIInt    pivot,t1,t2;
4324d615ea0SBarry Smith 
4334d615ea0SBarry Smith   PetscFunctionBegin;
434*59e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithArray,X,Y,n,pivot,t1,t2,ierr);
4355569a785SJunchao Zhang   PetscFunctionReturn(0);
4365569a785SJunchao Zhang }
4375569a785SJunchao Zhang 
4385569a785SJunchao Zhang /*@
4395569a785SJunchao Zhang    PetscSortMPIIntWithIntArray - Sorts an array of MPI integers in place in increasing order;
4405569a785SJunchao Zhang        changes a second array of Petsc intergers to match the sorted first array.
4415569a785SJunchao Zhang 
4425569a785SJunchao Zhang    Not Collective
4435569a785SJunchao Zhang 
4445569a785SJunchao Zhang    Input Parameters:
4455569a785SJunchao Zhang +  n  - number of values
446*59e16d97SJunchao Zhang .  X  - array of MPI integers
447*59e16d97SJunchao Zhang -  Y  - second array of Petsc integers
4485569a785SJunchao Zhang 
4495569a785SJunchao Zhang    Level: intermediate
4505569a785SJunchao Zhang 
4515569a785SJunchao Zhang    Notes: this routine is useful when one needs to sort MPI ranks with other integer arrays.
4525569a785SJunchao Zhang 
4535569a785SJunchao Zhang .seealso: PetscSortMPIIntWithArray()
4545569a785SJunchao Zhang @*/
455*59e16d97SJunchao Zhang PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n,PetscMPIInt *X,PetscInt *Y)
4565569a785SJunchao Zhang {
4575569a785SJunchao Zhang   PetscErrorCode ierr;
458*59e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
4595569a785SJunchao Zhang   PetscInt       t2;
4605569a785SJunchao Zhang 
4615569a785SJunchao Zhang   PetscFunctionBegin;
462*59e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithIntArray,X,Y,n,pivot,t1,t2,ierr);
463e5c89e4eSSatish Balay   PetscFunctionReturn(0);
464e5c89e4eSSatish Balay }
465e5c89e4eSSatish Balay 
466e5c89e4eSSatish Balay /*@
467e5c89e4eSSatish Balay    PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
468e5c89e4eSSatish Balay        changes a second SCALAR array to match the sorted first INTEGER array.
469e5c89e4eSSatish Balay 
470e5c89e4eSSatish Balay    Not Collective
471e5c89e4eSSatish Balay 
472e5c89e4eSSatish Balay    Input Parameters:
473e5c89e4eSSatish Balay +  n  - number of values
474*59e16d97SJunchao Zhang .  X  - array of integers
475*59e16d97SJunchao Zhang -  Y  - second array of scalars
476e5c89e4eSSatish Balay 
477e5c89e4eSSatish Balay    Level: intermediate
478e5c89e4eSSatish Balay 
479e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
480e5c89e4eSSatish Balay @*/
481*59e16d97SJunchao Zhang PetscErrorCode  PetscSortIntWithScalarArray(PetscInt n,PetscInt *X,PetscScalar *Y)
482e5c89e4eSSatish Balay {
483e5c89e4eSSatish Balay   PetscErrorCode ierr;
484*59e16d97SJunchao Zhang   PetscInt       pivot,t1;
485*59e16d97SJunchao Zhang   PetscScalar    t2;
486e5c89e4eSSatish Balay 
487e5c89e4eSSatish Balay   PetscFunctionBegin;
488*59e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithScalarArray,X,Y,n,pivot,t1,t2,ierr);
48917df18f2SToby Isaac   PetscFunctionReturn(0);
49017df18f2SToby Isaac }
49117df18f2SToby Isaac 
4926c2863d0SBarry Smith /*@C
49317df18f2SToby Isaac    PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order;
49417df18f2SToby Isaac        changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
49517df18f2SToby Isaac        provide workspace (the size of an element in the data array) to use when sorting.
49617df18f2SToby Isaac 
49717df18f2SToby Isaac    Not Collective
49817df18f2SToby Isaac 
49917df18f2SToby Isaac    Input Parameters:
50017df18f2SToby Isaac +  n  - number of values
501*59e16d97SJunchao Zhang .  X  - array of integers
502*59e16d97SJunchao Zhang .  Y  - second array of data
50317df18f2SToby Isaac .  size - sizeof elements in the data array in bytes
504*59e16d97SJunchao Zhang -  t2   - workspace of "size" bytes used when sorting
50517df18f2SToby Isaac 
50617df18f2SToby Isaac    Level: intermediate
50717df18f2SToby Isaac 
50817df18f2SToby Isaac .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
50917df18f2SToby Isaac @*/
510*59e16d97SJunchao Zhang PetscErrorCode  PetscSortIntWithDataArray(PetscInt n,PetscInt *X,void *Y,size_t size,void *t2)
51117df18f2SToby Isaac {
51217df18f2SToby Isaac   PetscErrorCode ierr;
513*59e16d97SJunchao Zhang   char           *YY = (char*)Y;
514*59e16d97SJunchao Zhang   PetscInt       i,j,p,t1,pivot,hi=n-1,l,r;
51517df18f2SToby Isaac 
51617df18f2SToby Isaac   PetscFunctionBegin;
51717df18f2SToby Isaac   if (n<8) {
518*59e16d97SJunchao Zhang     for (i=0; i<n; i++) {
519*59e16d97SJunchao Zhang       pivot = X[i];
520*59e16d97SJunchao Zhang       for (j=i+1; j<n; j++) {
521*59e16d97SJunchao Zhang         if (pivot > X[j]) {
522*59e16d97SJunchao Zhang           SWAP2Data(X[i],X[j],YY+size*i,YY+size*j,t1,t2,size);
523*59e16d97SJunchao Zhang           pivot = X[i];
52417df18f2SToby Isaac         }
52517df18f2SToby Isaac       }
52617df18f2SToby Isaac     }
52717df18f2SToby Isaac   } else {
528*59e16d97SJunchao Zhang     /* Two way partition */
529*59e16d97SJunchao Zhang     p     = MEDIAN(X,hi);
530*59e16d97SJunchao Zhang     pivot = X[p];
531*59e16d97SJunchao Zhang     l     = 0;
532*59e16d97SJunchao Zhang     r     = hi;
533*59e16d97SJunchao Zhang     while(1) {
534*59e16d97SJunchao Zhang       while (X[l] < pivot) l++;
535*59e16d97SJunchao Zhang       while (X[r] > pivot) r--;
536*59e16d97SJunchao Zhang       if (l >= r) {r++; break;}
537*59e16d97SJunchao Zhang       SWAP2Data(X[l],X[r],YY+size*l,YY+size*r,t1,t2,size);
538*59e16d97SJunchao Zhang       l++;
539*59e16d97SJunchao Zhang       r--;
540*59e16d97SJunchao Zhang     }
541*59e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(l,X,Y,size,t2);CHKERRQ(ierr);
542*59e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(hi-r+1,X+r,YY+size*r,size,t2);CHKERRQ(ierr);
54317df18f2SToby Isaac   }
54417df18f2SToby Isaac   PetscFunctionReturn(0);
54517df18f2SToby Isaac }
54617df18f2SToby Isaac 
54721e72a00SBarry Smith /*@
54821e72a00SBarry Smith    PetscMergeIntArray -     Merges two SORTED integer arrays, removes duplicate elements.
54921e72a00SBarry Smith 
55021e72a00SBarry Smith    Not Collective
55121e72a00SBarry Smith 
55221e72a00SBarry Smith    Input Parameters:
55321e72a00SBarry Smith +  an  - number of values in the first array
55421e72a00SBarry Smith .  aI  - first sorted array of integers
55521e72a00SBarry Smith .  bn  - number of values in the second array
55621e72a00SBarry Smith -  bI  - second array of integers
55721e72a00SBarry Smith 
55821e72a00SBarry Smith    Output Parameters:
55921e72a00SBarry Smith +  n   - number of values in the merged array
5606c2863d0SBarry Smith -  L   - merged sorted array, this is allocated if an array is not provided
56121e72a00SBarry Smith 
56221e72a00SBarry Smith    Level: intermediate
56321e72a00SBarry Smith 
56421e72a00SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
56521e72a00SBarry Smith @*/
5666c2863d0SBarry Smith PetscErrorCode  PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
56721e72a00SBarry Smith {
56821e72a00SBarry Smith   PetscErrorCode ierr;
56921e72a00SBarry Smith   PetscInt       *L_ = *L, ak, bk, k;
57021e72a00SBarry Smith 
57121e72a00SBarry Smith   if (!L_) {
57221e72a00SBarry Smith     ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr);
57321e72a00SBarry Smith     L_   = *L;
57421e72a00SBarry Smith   }
57521e72a00SBarry Smith   k = ak = bk = 0;
57621e72a00SBarry Smith   while (ak < an && bk < bn) {
57721e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
57821e72a00SBarry Smith       L_[k] = aI[ak];
57921e72a00SBarry Smith       ++ak;
58021e72a00SBarry Smith       ++bk;
58121e72a00SBarry Smith       ++k;
58221e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
58321e72a00SBarry Smith       L_[k] = aI[ak];
58421e72a00SBarry Smith       ++ak;
58521e72a00SBarry Smith       ++k;
58621e72a00SBarry Smith     } else {
58721e72a00SBarry Smith       L_[k] = bI[bk];
58821e72a00SBarry Smith       ++bk;
58921e72a00SBarry Smith       ++k;
59021e72a00SBarry Smith     }
59121e72a00SBarry Smith   }
59221e72a00SBarry Smith   if (ak < an) {
593580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
59421e72a00SBarry Smith     k   += (an-ak);
59521e72a00SBarry Smith   }
59621e72a00SBarry Smith   if (bk < bn) {
597580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
59821e72a00SBarry Smith     k   += (bn-bk);
59921e72a00SBarry Smith   }
60021e72a00SBarry Smith   *n = k;
60121e72a00SBarry Smith   PetscFunctionReturn(0);
60221e72a00SBarry Smith }
603b4301105SBarry Smith 
604c1f0200aSDmitry Karpeev /*@
60521e72a00SBarry Smith    PetscMergeIntArrayPair -     Merges two SORTED integer arrays that share NO common values along with an additional array of integers.
606c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
607c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
608c1f0200aSDmitry Karpeev 
609c1f0200aSDmitry Karpeev    Not Collective
610c1f0200aSDmitry Karpeev 
611c1f0200aSDmitry Karpeev    Input Parameters:
612c1f0200aSDmitry Karpeev +  an  - number of values in the first array
613c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
614c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
615c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
616c1f0200aSDmitry Karpeev .  bI  - second array of integers
617c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
618c1f0200aSDmitry Karpeev 
619c1f0200aSDmitry Karpeev    Output Parameters:
620c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
62114c49006SJed Brown .  L   - merged sorted array
622c1f0200aSDmitry Karpeev -  J   - merged additional array
623c1f0200aSDmitry Karpeev 
62495452b02SPatrick Sanan    Notes:
62595452b02SPatrick Sanan     if L or J point to non-null arrays then this routine will assume they are of the approproate size and use them, otherwise this routine will allocate space for them
626c1f0200aSDmitry Karpeev    Level: intermediate
627c1f0200aSDmitry Karpeev 
628c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
629c1f0200aSDmitry Karpeev @*/
6306c2863d0SBarry Smith PetscErrorCode  PetscMergeIntArrayPair(PetscInt an,const PetscInt aI[], const PetscInt aJ[], PetscInt bn, const PetscInt bI[], const PetscInt bJ[], PetscInt *n, PetscInt **L, PetscInt **J)
631c1f0200aSDmitry Karpeev {
632c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
63370d8d27cSBarry Smith   PetscInt       n_, *L_, *J_, ak, bk, k;
634c1f0200aSDmitry Karpeev 
63514c49006SJed Brown   PetscFunctionBegin;
63670d8d27cSBarry Smith   PetscValidIntPointer(L,8);
63770d8d27cSBarry Smith   PetscValidIntPointer(J,9);
638c1f0200aSDmitry Karpeev   n_ = an + bn;
639c1f0200aSDmitry Karpeev   *n = n_;
64070d8d27cSBarry Smith   if (!*L) {
641785e854fSJed Brown     ierr = PetscMalloc1(n_, L);CHKERRQ(ierr);
64270d8d27cSBarry Smith   }
643d7aa01a8SHong Zhang   L_ = *L;
64470d8d27cSBarry Smith   if (!*J) {
64570d8d27cSBarry Smith     ierr = PetscMalloc1(n_, J);CHKERRQ(ierr);
646c1f0200aSDmitry Karpeev   }
647c1f0200aSDmitry Karpeev   J_   = *J;
648c1f0200aSDmitry Karpeev   k = ak = bk = 0;
649c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
650c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
651d7aa01a8SHong Zhang       L_[k] = aI[ak];
652c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
653c1f0200aSDmitry Karpeev       ++ak;
654c1f0200aSDmitry Karpeev       ++k;
655a297a907SKarl Rupp     } else {
656d7aa01a8SHong Zhang       L_[k] = bI[bk];
657c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
658c1f0200aSDmitry Karpeev       ++bk;
659c1f0200aSDmitry Karpeev       ++k;
660c1f0200aSDmitry Karpeev     }
661c1f0200aSDmitry Karpeev   }
662c1f0200aSDmitry Karpeev   if (ak < an) {
663580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
664580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,aJ+ak,an-ak);CHKERRQ(ierr);
665c1f0200aSDmitry Karpeev     k   += (an-ak);
666c1f0200aSDmitry Karpeev   }
667c1f0200aSDmitry Karpeev   if (bk < bn) {
668580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
669580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,bJ+bk,bn-bk);CHKERRQ(ierr);
670c1f0200aSDmitry Karpeev   }
671c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
672c1f0200aSDmitry Karpeev }
673c1f0200aSDmitry Karpeev 
674e498c390SJed Brown /*@
675e498c390SJed Brown    PetscMergeMPIIntArray -     Merges two SORTED integer arrays.
676e498c390SJed Brown 
677e498c390SJed Brown    Not Collective
678e498c390SJed Brown 
679e498c390SJed Brown    Input Parameters:
680e498c390SJed Brown +  an  - number of values in the first array
681e498c390SJed Brown .  aI  - first sorted array of integers
682e498c390SJed Brown .  bn  - number of values in the second array
683e498c390SJed Brown -  bI  - second array of integers
684e498c390SJed Brown 
685e498c390SJed Brown    Output Parameters:
686e498c390SJed Brown +  n   - number of values in the merged array (<= an + bn)
687e498c390SJed Brown -  L   - merged sorted array, allocated if address of NULL pointer is passed
688e498c390SJed Brown 
689e498c390SJed Brown    Level: intermediate
690e498c390SJed Brown 
691e498c390SJed Brown .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
692e498c390SJed Brown @*/
693e498c390SJed Brown PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L)
694e498c390SJed Brown {
695e498c390SJed Brown   PetscErrorCode ierr;
696e498c390SJed Brown   PetscInt       ai,bi,k;
697e498c390SJed Brown 
698e498c390SJed Brown   PetscFunctionBegin;
699e498c390SJed Brown   if (!*L) {ierr = PetscMalloc1((an+bn),L);CHKERRQ(ierr);}
700e498c390SJed Brown   for (ai=0,bi=0,k=0; ai<an || bi<bn; ) {
701e498c390SJed Brown     PetscInt t = -1;
702e498c390SJed Brown     for ( ; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
703e498c390SJed Brown     for ( ; bi<bn && bI[bi] == t; bi++);
704e498c390SJed Brown     for ( ; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
705e498c390SJed Brown     for ( ; ai<an && aI[ai] == t; ai++);
706e498c390SJed Brown   }
707e498c390SJed Brown   *n = k;
708e498c390SJed Brown   PetscFunctionReturn(0);
709e498c390SJed Brown }
710e5c89e4eSSatish Balay 
7116c2863d0SBarry Smith /*@C
71242eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
71342eaa7f3SBarry Smith 
71442eaa7f3SBarry Smith    Not Collective
71542eaa7f3SBarry Smith 
71642eaa7f3SBarry Smith    Input Parameters:
71742eaa7f3SBarry Smith +  n  - number of values
71842eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
71942eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
72042eaa7f3SBarry Smith 
72142eaa7f3SBarry Smith    Output Parameters:
72242eaa7f3SBarry Smith +  Nlevels - the number of levels
72342eaa7f3SBarry Smith .  Level - for each node tells its level
72442eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
72542eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
72642eaa7f3SBarry Smith -  Column - for each id tells its column index
72742eaa7f3SBarry Smith 
7286c2863d0SBarry Smith    Level: developer
72942eaa7f3SBarry Smith 
73095452b02SPatrick Sanan    Notes:
73195452b02SPatrick Sanan     This code is not currently used
73242eaa7f3SBarry Smith 
73342eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation()
73442eaa7f3SBarry Smith @*/
7357087cfbeSBarry Smith PetscErrorCode  PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column)
73642eaa7f3SBarry Smith {
73742eaa7f3SBarry Smith   PetscInt       i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column;
73842eaa7f3SBarry Smith   PetscErrorCode ierr;
739ace3abfcSBarry Smith   PetscBool      done = PETSC_FALSE;
74042eaa7f3SBarry Smith 
74142eaa7f3SBarry Smith   PetscFunctionBegin;
74242eaa7f3SBarry Smith   if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set");
74342eaa7f3SBarry Smith   for (i=0; i<n; i++) {
74442eaa7f3SBarry Smith     if (mask[i]) continue;
74542eaa7f3SBarry Smith     if (parentid[i]  == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent");
74642eaa7f3SBarry Smith     if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked");
74742eaa7f3SBarry Smith   }
74842eaa7f3SBarry Smith 
74942eaa7f3SBarry Smith   for (i=0; i<n; i++) {
75042eaa7f3SBarry Smith     if (!mask[i]) nmask++;
75142eaa7f3SBarry Smith   }
75242eaa7f3SBarry Smith 
75342eaa7f3SBarry Smith   /* determine the level in the tree of each node */
7541795a4d1SJed Brown   ierr = PetscCalloc1(n,&level);CHKERRQ(ierr);
755a297a907SKarl Rupp 
75642eaa7f3SBarry Smith   level[0] = 1;
75742eaa7f3SBarry Smith   while (!done) {
75842eaa7f3SBarry Smith     done = PETSC_TRUE;
75942eaa7f3SBarry Smith     for (i=0; i<n; i++) {
76042eaa7f3SBarry Smith       if (mask[i]) continue;
76142eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
76242eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
76342eaa7f3SBarry Smith     }
76442eaa7f3SBarry Smith   }
76542eaa7f3SBarry Smith   for (i=0; i<n; i++) {
76642eaa7f3SBarry Smith     level[i]--;
76742eaa7f3SBarry Smith     nlevels = PetscMax(nlevels,level[i]);
76842eaa7f3SBarry Smith   }
76942eaa7f3SBarry Smith 
77042eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
7711795a4d1SJed Brown   ierr = PetscCalloc1(nlevels,&levelcnt);CHKERRQ(ierr);
77242eaa7f3SBarry Smith   for (i=0; i<n; i++) {
77342eaa7f3SBarry Smith     if (mask[i]) continue;
77442eaa7f3SBarry Smith     levelcnt[level[i]-1]++;
77542eaa7f3SBarry Smith   }
776a297a907SKarl Rupp   for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]);
77742eaa7f3SBarry Smith 
77842eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
779dcca6d9dSJed Brown   ierr = PetscMalloc2(levelmax,&workid,levelmax,&workparentid);CHKERRQ(ierr);
780785e854fSJed Brown   ierr = PetscMalloc1(nmask,&idbylevel);CHKERRQ(ierr);
78142eaa7f3SBarry Smith   for (j=1; j<=nlevels;j++) {
78242eaa7f3SBarry Smith     cnt = 0;
78342eaa7f3SBarry Smith     for (i=0; i<n; i++) {
78442eaa7f3SBarry Smith       if (mask[i]) continue;
78542eaa7f3SBarry Smith       if (level[i] != j) continue;
78642eaa7f3SBarry Smith       workid[cnt]         = i;
78742eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
78842eaa7f3SBarry Smith     }
78942eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
79042eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
79142eaa7f3SBarry Smith     ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr);
79242eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
79342eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
794580bdb30SBarry Smith     ierr  = PetscArraycpy(idbylevel+tcnt,workid,cnt);CHKERRQ(ierr);
79542eaa7f3SBarry Smith     tcnt += cnt;
79642eaa7f3SBarry Smith   }
79742eaa7f3SBarry Smith   if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes");
79842eaa7f3SBarry Smith   ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr);
79942eaa7f3SBarry Smith 
80042eaa7f3SBarry Smith   /* for each node list its column */
801785e854fSJed Brown   ierr = PetscMalloc1(n,&column);CHKERRQ(ierr);
80242eaa7f3SBarry Smith   cnt = 0;
80342eaa7f3SBarry Smith   for (j=0; j<nlevels; j++) {
80442eaa7f3SBarry Smith     for (i=0; i<levelcnt[j]; i++) {
80542eaa7f3SBarry Smith       column[idbylevel[cnt++]] = i;
80642eaa7f3SBarry Smith     }
80742eaa7f3SBarry Smith   }
80842eaa7f3SBarry Smith 
80942eaa7f3SBarry Smith   *Nlevels   = nlevels;
81042eaa7f3SBarry Smith   *Level     = level;
81142eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
81242eaa7f3SBarry Smith   *Idbylevel = idbylevel;
81342eaa7f3SBarry Smith   *Column    = column;
81442eaa7f3SBarry Smith   PetscFunctionReturn(0);
81542eaa7f3SBarry Smith }
816