xref: /petsc/src/sys/utils/sorti.c (revision 2a1da5286d8129285e2a5938dd041e517ba4f96c)
17d0a6c19SBarry Smith 
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file contains routines for sorting integers. Values are sorted in place.
459e16d97SJunchao 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*/
7f1cab4e1SJunchao Zhang #include <petsc/private/hashseti.h>
8e5c89e4eSSatish Balay 
9a8582168SJed Brown #define MEDIAN3(v,a,b,c)                                                        \
10a8582168SJed Brown   (v[a]<v[b]                                                                    \
11a8582168SJed Brown    ? (v[b]<v[c]                                                                 \
1259e16d97SJunchao Zhang       ? (b)                                                                     \
1359e16d97SJunchao Zhang       : (v[a]<v[c] ? (c) : (a)))                                                \
14a8582168SJed Brown    : (v[c]<v[b]                                                                 \
1559e16d97SJunchao Zhang       ? (b)                                                                     \
1659e16d97SJunchao Zhang       : (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 
2059e16d97SJunchao Zhang /* Swap one, two or three pairs. Each pair can have its own type */
2159e16d97SJunchao Zhang #define SWAP1(a,b,t1)               do {t1=a;a=b;b=t1;} while(0)
2259e16d97SJunchao Zhang #define SWAP2(a,b,c,d,t1,t2)        do {t1=a;a=b;b=t1; t2=c;c=d;d=t2;} while(0)
2359e16d97SJunchao 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)
2459e16d97SJunchao Zhang 
2559e16d97SJunchao Zhang /* Swap a & b, *c & *d. c, d, t2 are pointers to a type of size <siz> */
2659e16d97SJunchao Zhang #define SWAP2Data(a,b,c,d,t1,t2,siz)                                             \
2759e16d97SJunchao Zhang   do {                                                                           \
2859e16d97SJunchao Zhang     PetscErrorCode ierr;                                                         \
2959e16d97SJunchao Zhang     t1=a; a=b; b=t1;                                                             \
3059e16d97SJunchao Zhang     ierr = PetscMemcpy(t2,c,siz);CHKERRQ(ierr);                                  \
3159e16d97SJunchao Zhang     ierr = PetscMemcpy(c,d,siz);CHKERRQ(ierr);                                   \
3259e16d97SJunchao Zhang     ierr = PetscMemcpy(d,t2,siz);CHKERRQ(ierr);                                  \
3359e16d97SJunchao Zhang   } while(0)
34e5c89e4eSSatish Balay 
35e5c89e4eSSatish Balay /*
3659e16d97SJunchao Zhang    Partition X[lo,hi] into two parts: X[lo,l) <= pivot; X[r,hi] > pivot
37e5c89e4eSSatish Balay 
3859e16d97SJunchao Zhang    Input Parameters:
3959e16d97SJunchao Zhang     + X         - array to partition
4059e16d97SJunchao Zhang     . pivot     - a pivot of X[]
4159e16d97SJunchao Zhang     . t1        - temp variable for X
4259e16d97SJunchao Zhang     - lo,hi     - lower and upper bound of the array
4359e16d97SJunchao Zhang 
4459e16d97SJunchao Zhang    Output Parameters:
4559e16d97SJunchao Zhang     + l,r       - of type PetscInt
4659e16d97SJunchao Zhang 
4759e16d97SJunchao Zhang    Notes:
4859e16d97SJunchao Zhang     The TwoWayPartition2/3 variants also partition other arrays along with X.
4959e16d97SJunchao Zhang     These arrays can have different types, so they provide their own temp t2,t3
5059e16d97SJunchao Zhang  */
5159e16d97SJunchao Zhang #define TwoWayPartition1(X,pivot,t1,lo,hi,l,r)                                   \
5259e16d97SJunchao Zhang   do {                                                                           \
5359e16d97SJunchao Zhang     l = lo;                                                                      \
5459e16d97SJunchao Zhang     r = hi;                                                                      \
5559e16d97SJunchao Zhang     while(1) {                                                                   \
5659e16d97SJunchao Zhang       while (X[l] < pivot) l++;                                                  \
5759e16d97SJunchao Zhang       while (X[r] > pivot) r--;                                                  \
5859e16d97SJunchao Zhang       if (l >= r) {r++; break;}                                                  \
5959e16d97SJunchao Zhang       SWAP1(X[l],X[r],t1);                                                       \
6059e16d97SJunchao Zhang       l++;                                                                       \
6159e16d97SJunchao Zhang       r--;                                                                       \
6259e16d97SJunchao Zhang     }                                                                            \
6359e16d97SJunchao Zhang   } while(0)
6459e16d97SJunchao Zhang 
65ce605777SToby Isaac /*
66ce605777SToby Isaac    Partition X[lo,hi] into two parts: X[lo,l) >= pivot; X[r,hi] < pivot
67ce605777SToby Isaac 
68ce605777SToby Isaac    Input Parameters:
69ce605777SToby Isaac     + X         - array to partition
70ce605777SToby Isaac     . pivot     - a pivot of X[]
71ce605777SToby Isaac     . t1        - temp variable for X
72ce605777SToby Isaac     - lo,hi     - lower and upper bound of the array
73ce605777SToby Isaac 
74ce605777SToby Isaac    Output Parameters:
75ce605777SToby Isaac     + l,r       - of type PetscInt
76ce605777SToby Isaac 
77ce605777SToby Isaac    Notes:
78ce605777SToby Isaac     The TwoWayPartition2/3 variants also partition other arrays along with X.
79ce605777SToby Isaac     These arrays can have different types, so they provide their own temp t2,t3
80ce605777SToby Isaac  */
81ce605777SToby Isaac #define TwoWayPartitionReverse1(X,pivot,t1,lo,hi,l,r)                                   \
82ce605777SToby Isaac   do {                                                                           \
83ce605777SToby Isaac     l = lo;                                                                      \
84ce605777SToby Isaac     r = hi;                                                                      \
85ce605777SToby Isaac     while(1) {                                                                   \
86ce605777SToby Isaac       while (X[l] > pivot) l++;                                                  \
87ce605777SToby Isaac       while (X[r] < pivot) r--;                                                  \
88ce605777SToby Isaac       if (l >= r) {r++; break;}                                                  \
89ce605777SToby Isaac       SWAP1(X[l],X[r],t1);                                                       \
90ce605777SToby Isaac       l++;                                                                       \
91ce605777SToby Isaac       r--;                                                                       \
92ce605777SToby Isaac     }                                                                            \
93ce605777SToby Isaac   } while(0)
94ce605777SToby Isaac 
9559e16d97SJunchao Zhang #define TwoWayPartition2(X,Y,pivot,t1,t2,lo,hi,l,r)                              \
9659e16d97SJunchao Zhang   do {                                                                           \
9759e16d97SJunchao Zhang     l = lo;                                                                      \
9859e16d97SJunchao Zhang     r = hi;                                                                      \
9959e16d97SJunchao Zhang     while(1) {                                                                   \
10059e16d97SJunchao Zhang       while (X[l] < pivot) l++;                                                  \
10159e16d97SJunchao Zhang       while (X[r] > pivot) r--;                                                  \
10259e16d97SJunchao Zhang       if (l >= r) {r++; break;}                                                  \
10359e16d97SJunchao Zhang       SWAP2(X[l],X[r],Y[l],Y[r],t1,t2);                                          \
10459e16d97SJunchao Zhang       l++;                                                                       \
10559e16d97SJunchao Zhang       r--;                                                                       \
10659e16d97SJunchao Zhang     }                                                                            \
10759e16d97SJunchao Zhang   } while(0)
10859e16d97SJunchao Zhang 
10959e16d97SJunchao Zhang #define TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,lo,hi,l,r)                         \
11059e16d97SJunchao Zhang   do {                                                                           \
11159e16d97SJunchao Zhang     l = lo;                                                                      \
11259e16d97SJunchao Zhang     r = hi;                                                                      \
11359e16d97SJunchao Zhang     while(1) {                                                                   \
11459e16d97SJunchao Zhang       while (X[l] < pivot) l++;                                                  \
11559e16d97SJunchao Zhang       while (X[r] > pivot) r--;                                                  \
11659e16d97SJunchao Zhang       if (l >= r) {r++; break;}                                                  \
11759e16d97SJunchao Zhang       SWAP3(X[l],X[r],Y[l],Y[r],Z[l],Z[r],t1,t2,t3);                             \
11859e16d97SJunchao Zhang       l++;                                                                       \
11959e16d97SJunchao Zhang       r--;                                                                       \
12059e16d97SJunchao Zhang     }                                                                            \
12159e16d97SJunchao Zhang   } while(0)
12259e16d97SJunchao Zhang 
12359e16d97SJunchao Zhang /* Templates for similar functions used below */
12459e16d97SJunchao Zhang #define QuickSort1(FuncName,X,n,pivot,t1,ierr)                                   \
12559e16d97SJunchao Zhang   do {                                                                           \
12659e16d97SJunchao Zhang     PetscInt i,j,p,l,r,hi=n-1;                                                   \
12759e16d97SJunchao Zhang     if (n < 8) {                                                                 \
12859e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
12959e16d97SJunchao Zhang         pivot = X[i];                                                            \
13059e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
13159e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
13259e16d97SJunchao Zhang             SWAP1(X[i],X[j],t1);                                                 \
13359e16d97SJunchao Zhang             pivot = X[i];                                                        \
13459e16d97SJunchao Zhang           }                                                                      \
13559e16d97SJunchao Zhang         }                                                                        \
13659e16d97SJunchao Zhang       }                                                                          \
13759e16d97SJunchao Zhang     } else {                                                                     \
13859e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
13959e16d97SJunchao Zhang       pivot = X[p];                                                              \
14059e16d97SJunchao Zhang       TwoWayPartition1(X,pivot,t1,0,hi,l,r);                                     \
14159e16d97SJunchao Zhang       ierr  = FuncName(l,X);CHKERRQ(ierr);                                       \
14259e16d97SJunchao Zhang       ierr  = FuncName(hi-r+1,X+r);CHKERRQ(ierr);                                \
14359e16d97SJunchao Zhang     }                                                                            \
14459e16d97SJunchao Zhang   } while(0)
14559e16d97SJunchao Zhang 
146ce605777SToby Isaac /* Templates for similar functions used below */
147ce605777SToby Isaac #define QuickSortReverse1(FuncName,X,n,pivot,t1,ierr)                            \
148ce605777SToby Isaac   do {                                                                           \
149ce605777SToby Isaac     PetscInt i,j,p,l,r,hi=n-1;                                                   \
150ce605777SToby Isaac     if (n < 8) {                                                                 \
151ce605777SToby Isaac       for (i=0; i<n; i++) {                                                      \
152ce605777SToby Isaac         pivot = X[i];                                                            \
153ce605777SToby Isaac         for (j=i+1; j<n; j++) {                                                  \
154ce605777SToby Isaac           if (pivot < X[j]) {                                                    \
155ce605777SToby Isaac             SWAP1(X[i],X[j],t1);                                                 \
156ce605777SToby Isaac             pivot = X[i];                                                        \
157ce605777SToby Isaac           }                                                                      \
158ce605777SToby Isaac         }                                                                        \
159ce605777SToby Isaac       }                                                                          \
160ce605777SToby Isaac     } else {                                                                     \
161ce605777SToby Isaac       p     = MEDIAN(X,hi);                                                      \
162ce605777SToby Isaac       pivot = X[p];                                                              \
163ce605777SToby Isaac       TwoWayPartitionReverse1(X,pivot,t1,0,hi,l,r);                              \
164ce605777SToby Isaac       ierr  = FuncName(l,X);CHKERRQ(ierr);                                       \
165ce605777SToby Isaac       ierr  = FuncName(hi-r+1,X+r);CHKERRQ(ierr);                                \
166ce605777SToby Isaac     }                                                                            \
167ce605777SToby Isaac   } while(0)
168ce605777SToby Isaac 
16959e16d97SJunchao Zhang #define QuickSort2(FuncName,X,Y,n,pivot,t1,t2,ierr)                              \
17059e16d97SJunchao Zhang   do {                                                                           \
17159e16d97SJunchao Zhang     PetscInt i,j,p,l,r,hi=n-1;                                                   \
17259e16d97SJunchao Zhang     if (n < 8) {                                                                 \
17359e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
17459e16d97SJunchao Zhang         pivot = X[i];                                                            \
17559e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
17659e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
17759e16d97SJunchao Zhang             SWAP2(X[i],X[j],Y[i],Y[j],t1,t2);                                    \
17859e16d97SJunchao Zhang             pivot = X[i];                                                        \
17959e16d97SJunchao Zhang           }                                                                      \
18059e16d97SJunchao Zhang         }                                                                        \
18159e16d97SJunchao Zhang       }                                                                          \
18259e16d97SJunchao Zhang     } else {                                                                     \
18359e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
18459e16d97SJunchao Zhang       pivot = X[p];                                                              \
18559e16d97SJunchao Zhang       TwoWayPartition2(X,Y,pivot,t1,t2,0,hi,l,r);                                \
18659e16d97SJunchao Zhang       ierr  = FuncName(l,X,Y);CHKERRQ(ierr);                                     \
18759e16d97SJunchao Zhang       ierr  = FuncName(hi-r+1,X+r,Y+r);CHKERRQ(ierr);                            \
18859e16d97SJunchao Zhang     }                                                                            \
18959e16d97SJunchao Zhang   } while(0)
19059e16d97SJunchao Zhang 
19159e16d97SJunchao Zhang #define QuickSort3(FuncName,X,Y,Z,n,pivot,t1,t2,t3,ierr)                         \
19259e16d97SJunchao Zhang   do {                                                                           \
19359e16d97SJunchao Zhang     PetscInt i,j,p,l,r,hi=n-1;                                                   \
19459e16d97SJunchao Zhang     if (n < 8) {                                                                 \
19559e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
19659e16d97SJunchao Zhang         pivot = X[i];                                                            \
19759e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
19859e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
19959e16d97SJunchao Zhang             SWAP3(X[i],X[j],Y[i],Y[j],Z[i],Z[j],t1,t2,t3);                       \
20059e16d97SJunchao Zhang             pivot = X[i];                                                        \
20159e16d97SJunchao Zhang           }                                                                      \
20259e16d97SJunchao Zhang         }                                                                        \
20359e16d97SJunchao Zhang       }                                                                          \
20459e16d97SJunchao Zhang     } else {                                                                     \
20559e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
20659e16d97SJunchao Zhang       pivot = X[p];                                                              \
20759e16d97SJunchao Zhang       TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,0,hi,l,r);                           \
20859e16d97SJunchao Zhang       ierr  = FuncName(l,X,Y,Z);CHKERRQ(ierr);                                   \
20959e16d97SJunchao Zhang       ierr  = FuncName(hi-r+1,X+r,Y+r,Z+r);CHKERRQ(ierr);                        \
21059e16d97SJunchao Zhang     }                                                                            \
21159e16d97SJunchao Zhang   } while(0)
212e5c89e4eSSatish Balay 
213e5c89e4eSSatish Balay /*@
214e5c89e4eSSatish Balay    PetscSortInt - Sorts an array of integers in place in increasing order.
215e5c89e4eSSatish Balay 
216e5c89e4eSSatish Balay    Not Collective
217e5c89e4eSSatish Balay 
218e5c89e4eSSatish Balay    Input Parameters:
219e5c89e4eSSatish Balay +  n  - number of values
22059e16d97SJunchao Zhang -  X  - array of integers
221e5c89e4eSSatish Balay 
222e5c89e4eSSatish Balay    Level: intermediate
223e5c89e4eSSatish Balay 
224e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntWithPermutation()
225e5c89e4eSSatish Balay @*/
226b4f531b9SJunchao Zhang PetscErrorCode  PetscSortInt(PetscInt n,PetscInt X[])
227e5c89e4eSSatish Balay {
22859e16d97SJunchao Zhang   PetscErrorCode ierr;
22959e16d97SJunchao Zhang   PetscInt       pivot,t1;
230e5c89e4eSSatish Balay 
231e5c89e4eSSatish Balay   PetscFunctionBegin;
23259e16d97SJunchao Zhang   QuickSort1(PetscSortInt,X,n,pivot,t1,ierr);
233e5c89e4eSSatish Balay   PetscFunctionReturn(0);
234e5c89e4eSSatish Balay }
235e5c89e4eSSatish Balay 
2361db36a52SBarry Smith /*@
237ce605777SToby Isaac    PetscSortReverseInt - Sorts an array of integers in place in decreasing order.
238ce605777SToby Isaac 
239ce605777SToby Isaac    Not Collective
240ce605777SToby Isaac 
241ce605777SToby Isaac    Input Parameters:
242ce605777SToby Isaac +  n  - number of values
243ce605777SToby Isaac -  X  - array of integers
244ce605777SToby Isaac 
245ce605777SToby Isaac    Level: intermediate
246ce605777SToby Isaac 
247ce605777SToby Isaac .seealso: PetscSortInt(), PetscSortIntWithPermutation()
248ce605777SToby Isaac @*/
249ce605777SToby Isaac PetscErrorCode  PetscSortReverseInt(PetscInt n,PetscInt X[])
250ce605777SToby Isaac {
251ce605777SToby Isaac   PetscErrorCode ierr;
252ce605777SToby Isaac   PetscInt       pivot,t1;
253ce605777SToby Isaac 
254ce605777SToby Isaac   PetscFunctionBegin;
255ce605777SToby Isaac   QuickSortReverse1(PetscSortReverseInt,X,n,pivot,t1,ierr);
256ce605777SToby Isaac   PetscFunctionReturn(0);
257ce605777SToby Isaac }
258ce605777SToby Isaac 
259ce605777SToby Isaac /*@
26022ab5688SLisandro Dalcin    PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted input array
26122ab5688SLisandro Dalcin 
26222ab5688SLisandro Dalcin    Not Collective
26322ab5688SLisandro Dalcin 
26422ab5688SLisandro Dalcin    Input Parameters:
26522ab5688SLisandro Dalcin +  n  - number of values
26659e16d97SJunchao Zhang -  X  - sorted array of integers
26722ab5688SLisandro Dalcin 
26822ab5688SLisandro Dalcin    Output Parameter:
26922ab5688SLisandro Dalcin .  n - number of non-redundant values
27022ab5688SLisandro Dalcin 
27122ab5688SLisandro Dalcin    Level: intermediate
27222ab5688SLisandro Dalcin 
27322ab5688SLisandro Dalcin .seealso: PetscSortInt()
27422ab5688SLisandro Dalcin @*/
27559e16d97SJunchao Zhang PetscErrorCode  PetscSortedRemoveDupsInt(PetscInt *n,PetscInt X[])
27622ab5688SLisandro Dalcin {
27722ab5688SLisandro Dalcin   PetscInt i,s = 0,N = *n, b = 0;
27822ab5688SLisandro Dalcin 
27922ab5688SLisandro Dalcin   PetscFunctionBegin;
28022ab5688SLisandro Dalcin   for (i=0; i<N-1; i++) {
28159e16d97SJunchao Zhang     if (PetscUnlikely(X[b+s+1] < X[b])) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Input array is not sorted");
28259e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
28359e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
28422ab5688SLisandro Dalcin     } else s++;
28522ab5688SLisandro Dalcin   }
28622ab5688SLisandro Dalcin   *n = N - s;
28722ab5688SLisandro Dalcin   PetscFunctionReturn(0);
28822ab5688SLisandro Dalcin }
28922ab5688SLisandro Dalcin 
29022ab5688SLisandro Dalcin /*@
2911db36a52SBarry Smith    PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries
2921db36a52SBarry Smith 
2931db36a52SBarry Smith    Not Collective
2941db36a52SBarry Smith 
2951db36a52SBarry Smith    Input Parameters:
2961db36a52SBarry Smith +  n  - number of values
29759e16d97SJunchao Zhang -  X  - array of integers
2981db36a52SBarry Smith 
2991db36a52SBarry Smith    Output Parameter:
3001db36a52SBarry Smith .  n - number of non-redundant values
3011db36a52SBarry Smith 
3021db36a52SBarry Smith    Level: intermediate
3031db36a52SBarry Smith 
30422ab5688SLisandro Dalcin .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortedRemoveDupsInt()
3051db36a52SBarry Smith @*/
30659e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsInt(PetscInt *n,PetscInt X[])
3071db36a52SBarry Smith {
3081db36a52SBarry Smith   PetscErrorCode ierr;
3091db36a52SBarry Smith 
3101db36a52SBarry Smith   PetscFunctionBegin;
31159e16d97SJunchao Zhang   ierr = PetscSortInt(*n,X);CHKERRQ(ierr);
31259e16d97SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(n,X);CHKERRQ(ierr);
3131db36a52SBarry Smith   PetscFunctionReturn(0);
3141db36a52SBarry Smith }
3151db36a52SBarry Smith 
31660e03357SMatthew G Knepley /*@
317d9f1162dSJed Brown   PetscFindInt - Finds integer in a sorted array of integers
31860e03357SMatthew G Knepley 
31960e03357SMatthew G Knepley    Not Collective
32060e03357SMatthew G Knepley 
32160e03357SMatthew G Knepley    Input Parameters:
32260e03357SMatthew G Knepley +  key - the integer to locate
323d9f1162dSJed Brown .  n   - number of values in the array
32459e16d97SJunchao Zhang -  X  - array of integers
32560e03357SMatthew G Knepley 
32660e03357SMatthew G Knepley    Output Parameter:
327d9f1162dSJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
32860e03357SMatthew G Knepley 
32960e03357SMatthew G Knepley    Level: intermediate
33060e03357SMatthew G Knepley 
33160e03357SMatthew G Knepley .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
33260e03357SMatthew G Knepley @*/
33359e16d97SJunchao Zhang PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
33460e03357SMatthew G Knepley {
335d9f1162dSJed Brown   PetscInt lo = 0,hi = n;
33660e03357SMatthew G Knepley 
33760e03357SMatthew G Knepley   PetscFunctionBegin;
33860e03357SMatthew G Knepley   PetscValidPointer(loc,4);
33998781d82SMatthew G Knepley   if (!n) {*loc = -1; PetscFunctionReturn(0);}
34059e16d97SJunchao Zhang   PetscValidPointer(X,3);
341d9f1162dSJed Brown   while (hi - lo > 1) {
342d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo)/2;
34359e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
344d9f1162dSJed Brown     else               lo = mid;
34560e03357SMatthew G Knepley   }
34659e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
34760e03357SMatthew G Knepley   PetscFunctionReturn(0);
34860e03357SMatthew G Knepley }
34960e03357SMatthew G Knepley 
350d2aeb606SJed Brown /*@
351f1cab4e1SJunchao Zhang   PetscCheckDupsInt - Checks if an integer array has duplicates
352f1cab4e1SJunchao Zhang 
353f1cab4e1SJunchao Zhang    Not Collective
354f1cab4e1SJunchao Zhang 
355f1cab4e1SJunchao Zhang    Input Parameters:
356f1cab4e1SJunchao Zhang +  n  - number of values in the array
357f1cab4e1SJunchao Zhang -  X  - array of integers
358f1cab4e1SJunchao Zhang 
359f1cab4e1SJunchao Zhang 
360f1cab4e1SJunchao Zhang    Output Parameter:
361f1cab4e1SJunchao Zhang .  dups - True if the array has dups, otherwise false
362f1cab4e1SJunchao Zhang 
363f1cab4e1SJunchao Zhang    Level: intermediate
364f1cab4e1SJunchao Zhang 
365f1cab4e1SJunchao Zhang .seealso: PetscSortRemoveDupsInt()
366f1cab4e1SJunchao Zhang @*/
367f1cab4e1SJunchao Zhang PetscErrorCode PetscCheckDupsInt(PetscInt n,const PetscInt X[],PetscBool *dups)
368f1cab4e1SJunchao Zhang {
369f1cab4e1SJunchao Zhang   PetscErrorCode ierr;
370f1cab4e1SJunchao Zhang   PetscInt       i;
371f1cab4e1SJunchao Zhang   PetscHSetI     ht;
372f1cab4e1SJunchao Zhang   PetscBool      missing;
373f1cab4e1SJunchao Zhang 
374f1cab4e1SJunchao Zhang   PetscFunctionBegin;
375f1cab4e1SJunchao Zhang   PetscValidPointer(dups,3);
376f1cab4e1SJunchao Zhang   *dups = PETSC_FALSE;
377f1cab4e1SJunchao Zhang   if (n > 1) {
378f1cab4e1SJunchao Zhang     ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
379f1cab4e1SJunchao Zhang     ierr = PetscHSetIResize(ht,n);CHKERRQ(ierr);
380f1cab4e1SJunchao Zhang     for (i=0; i<n; i++) {
381f1cab4e1SJunchao Zhang       ierr = PetscHSetIQueryAdd(ht,X[i],&missing);CHKERRQ(ierr);
382f1cab4e1SJunchao Zhang       if(!missing) {*dups = PETSC_TRUE; break;}
383f1cab4e1SJunchao Zhang     }
384f1cab4e1SJunchao Zhang     ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
385f1cab4e1SJunchao Zhang   }
386f1cab4e1SJunchao Zhang   PetscFunctionReturn(0);
387f1cab4e1SJunchao Zhang }
388f1cab4e1SJunchao Zhang 
389f1cab4e1SJunchao Zhang /*@
390d2aeb606SJed Brown   PetscFindMPIInt - Finds MPI integer in a sorted array of integers
391d2aeb606SJed Brown 
392d2aeb606SJed Brown    Not Collective
393d2aeb606SJed Brown 
394d2aeb606SJed Brown    Input Parameters:
395d2aeb606SJed Brown +  key - the integer to locate
396d2aeb606SJed Brown .  n   - number of values in the array
39759e16d97SJunchao Zhang -  X   - array of integers
398d2aeb606SJed Brown 
399d2aeb606SJed Brown    Output Parameter:
400d2aeb606SJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
401d2aeb606SJed Brown 
402d2aeb606SJed Brown    Level: intermediate
403d2aeb606SJed Brown 
404d2aeb606SJed Brown .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
405d2aeb606SJed Brown @*/
40659e16d97SJunchao Zhang PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
407d2aeb606SJed Brown {
408d2aeb606SJed Brown   PetscInt lo = 0,hi = n;
409d2aeb606SJed Brown 
410d2aeb606SJed Brown   PetscFunctionBegin;
411d2aeb606SJed Brown   PetscValidPointer(loc,4);
412d2aeb606SJed Brown   if (!n) {*loc = -1; PetscFunctionReturn(0);}
41359e16d97SJunchao Zhang   PetscValidPointer(X,3);
414d2aeb606SJed Brown   while (hi - lo > 1) {
415d2aeb606SJed Brown     PetscInt mid = lo + (hi - lo)/2;
41659e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
417d2aeb606SJed Brown     else               lo = mid;
418d2aeb606SJed Brown   }
41959e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
420e5c89e4eSSatish Balay   PetscFunctionReturn(0);
421e5c89e4eSSatish Balay }
422e5c89e4eSSatish Balay 
423e5c89e4eSSatish Balay /*@
424e5c89e4eSSatish Balay    PetscSortIntWithArray - Sorts an array of integers in place in increasing order;
425e5c89e4eSSatish Balay        changes a second array to match the sorted first array.
426e5c89e4eSSatish Balay 
427e5c89e4eSSatish Balay    Not Collective
428e5c89e4eSSatish Balay 
429e5c89e4eSSatish Balay    Input Parameters:
430e5c89e4eSSatish Balay +  n  - number of values
43159e16d97SJunchao Zhang .  X  - array of integers
43259e16d97SJunchao Zhang -  Y  - second array of integers
433e5c89e4eSSatish Balay 
434e5c89e4eSSatish Balay    Level: intermediate
435e5c89e4eSSatish Balay 
436e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
437e5c89e4eSSatish Balay @*/
438b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithArray(PetscInt n,PetscInt X[],PetscInt Y[])
439e5c89e4eSSatish Balay {
440e5c89e4eSSatish Balay   PetscErrorCode ierr;
44159e16d97SJunchao Zhang   PetscInt       pivot,t1,t2;
442e5c89e4eSSatish Balay 
443e5c89e4eSSatish Balay   PetscFunctionBegin;
44459e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithArray,X,Y,n,pivot,t1,t2,ierr);
445c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
446c1f0200aSDmitry Karpeev }
447c1f0200aSDmitry Karpeev 
448c1f0200aSDmitry Karpeev /*@
449c1f0200aSDmitry Karpeev    PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order;
450c1f0200aSDmitry Karpeev        changes a pair of integer arrays to match the sorted first array.
451c1f0200aSDmitry Karpeev 
452c1f0200aSDmitry Karpeev    Not Collective
453c1f0200aSDmitry Karpeev 
454c1f0200aSDmitry Karpeev    Input Parameters:
455c1f0200aSDmitry Karpeev +  n  - number of values
45659e16d97SJunchao Zhang .  X  - array of integers
45759e16d97SJunchao Zhang .  Y  - second array of integers (first array of the pair)
45859e16d97SJunchao Zhang -  Z  - third array of integers  (second array of the pair)
459c1f0200aSDmitry Karpeev 
460c1f0200aSDmitry Karpeev    Level: intermediate
461c1f0200aSDmitry Karpeev 
462c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray()
463c1f0200aSDmitry Karpeev @*/
464b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithArrayPair(PetscInt n,PetscInt X[],PetscInt Y[],PetscInt Z[])
465c1f0200aSDmitry Karpeev {
466c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
46759e16d97SJunchao Zhang   PetscInt       pivot,t1,t2,t3;
468c1f0200aSDmitry Karpeev 
469c1f0200aSDmitry Karpeev   PetscFunctionBegin;
47059e16d97SJunchao Zhang   QuickSort3(PetscSortIntWithArrayPair,X,Y,Z,n,pivot,t1,t2,t3,ierr);
471c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
472c1f0200aSDmitry Karpeev }
473c1f0200aSDmitry Karpeev 
47417d7d925SStefano Zampini /*@
47517d7d925SStefano Zampini    PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order.
47617d7d925SStefano Zampini 
47717d7d925SStefano Zampini    Not Collective
47817d7d925SStefano Zampini 
47917d7d925SStefano Zampini    Input Parameters:
48017d7d925SStefano Zampini +  n  - number of values
48159e16d97SJunchao Zhang -  X  - array of integers
48217d7d925SStefano Zampini 
48317d7d925SStefano Zampini    Level: intermediate
48417d7d925SStefano Zampini 
48517d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation()
48617d7d925SStefano Zampini @*/
487b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIInt(PetscInt n,PetscMPIInt X[])
48817d7d925SStefano Zampini {
48959e16d97SJunchao Zhang   PetscErrorCode ierr;
49059e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
49117d7d925SStefano Zampini 
49217d7d925SStefano Zampini   PetscFunctionBegin;
49359e16d97SJunchao Zhang   QuickSort1(PetscSortMPIInt,X,n,pivot,t1,ierr);
49417d7d925SStefano Zampini   PetscFunctionReturn(0);
49517d7d925SStefano Zampini }
49617d7d925SStefano Zampini 
49717d7d925SStefano Zampini /*@
49817d7d925SStefano Zampini    PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries
49917d7d925SStefano Zampini 
50017d7d925SStefano Zampini    Not Collective
50117d7d925SStefano Zampini 
50217d7d925SStefano Zampini    Input Parameters:
50317d7d925SStefano Zampini +  n  - number of values
50459e16d97SJunchao Zhang -  X  - array of integers
50517d7d925SStefano Zampini 
50617d7d925SStefano Zampini    Output Parameter:
50717d7d925SStefano Zampini .  n - number of non-redundant values
50817d7d925SStefano Zampini 
50917d7d925SStefano Zampini    Level: intermediate
51017d7d925SStefano Zampini 
51117d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
51217d7d925SStefano Zampini @*/
51359e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt X[])
51417d7d925SStefano Zampini {
51517d7d925SStefano Zampini   PetscErrorCode ierr;
51617d7d925SStefano Zampini   PetscInt       i,s = 0,N = *n, b = 0;
51717d7d925SStefano Zampini 
51817d7d925SStefano Zampini   PetscFunctionBegin;
51959e16d97SJunchao Zhang   ierr = PetscSortMPIInt(N,X);CHKERRQ(ierr);
52017d7d925SStefano Zampini   for (i=0; i<N-1; i++) {
52159e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
52259e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
523a297a907SKarl Rupp     } else s++;
52417d7d925SStefano Zampini   }
52517d7d925SStefano Zampini   *n = N - s;
52617d7d925SStefano Zampini   PetscFunctionReturn(0);
52717d7d925SStefano Zampini }
528c1f0200aSDmitry Karpeev 
5294d615ea0SBarry Smith /*@
5304d615ea0SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order;
5314d615ea0SBarry Smith        changes a second array to match the sorted first array.
5324d615ea0SBarry Smith 
5334d615ea0SBarry Smith    Not Collective
5344d615ea0SBarry Smith 
5354d615ea0SBarry Smith    Input Parameters:
5364d615ea0SBarry Smith +  n  - number of values
53759e16d97SJunchao Zhang .  X  - array of integers
53859e16d97SJunchao Zhang -  Y  - second array of integers
5394d615ea0SBarry Smith 
5404d615ea0SBarry Smith    Level: intermediate
5414d615ea0SBarry Smith 
5424d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
5434d615ea0SBarry Smith @*/
544b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt X[],PetscMPIInt Y[])
5454d615ea0SBarry Smith {
5464d615ea0SBarry Smith   PetscErrorCode ierr;
54759e16d97SJunchao Zhang   PetscMPIInt    pivot,t1,t2;
5484d615ea0SBarry Smith 
5494d615ea0SBarry Smith   PetscFunctionBegin;
55059e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithArray,X,Y,n,pivot,t1,t2,ierr);
5515569a785SJunchao Zhang   PetscFunctionReturn(0);
5525569a785SJunchao Zhang }
5535569a785SJunchao Zhang 
5545569a785SJunchao Zhang /*@
5555569a785SJunchao Zhang    PetscSortMPIIntWithIntArray - Sorts an array of MPI integers in place in increasing order;
5565569a785SJunchao Zhang        changes a second array of Petsc intergers to match the sorted first array.
5575569a785SJunchao Zhang 
5585569a785SJunchao Zhang    Not Collective
5595569a785SJunchao Zhang 
5605569a785SJunchao Zhang    Input Parameters:
5615569a785SJunchao Zhang +  n  - number of values
56259e16d97SJunchao Zhang .  X  - array of MPI integers
56359e16d97SJunchao Zhang -  Y  - second array of Petsc integers
5645569a785SJunchao Zhang 
5655569a785SJunchao Zhang    Level: intermediate
5665569a785SJunchao Zhang 
5675569a785SJunchao Zhang    Notes: this routine is useful when one needs to sort MPI ranks with other integer arrays.
5685569a785SJunchao Zhang 
5695569a785SJunchao Zhang .seealso: PetscSortMPIIntWithArray()
5705569a785SJunchao Zhang @*/
571b4f531b9SJunchao Zhang PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n,PetscMPIInt X[],PetscInt Y[])
5725569a785SJunchao Zhang {
5735569a785SJunchao Zhang   PetscErrorCode ierr;
57459e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
5755569a785SJunchao Zhang   PetscInt       t2;
5765569a785SJunchao Zhang 
5775569a785SJunchao Zhang   PetscFunctionBegin;
57859e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithIntArray,X,Y,n,pivot,t1,t2,ierr);
579e5c89e4eSSatish Balay   PetscFunctionReturn(0);
580e5c89e4eSSatish Balay }
581e5c89e4eSSatish Balay 
582e5c89e4eSSatish Balay /*@
583e5c89e4eSSatish Balay    PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
584e5c89e4eSSatish Balay        changes a second SCALAR array to match the sorted first INTEGER array.
585e5c89e4eSSatish Balay 
586e5c89e4eSSatish Balay    Not Collective
587e5c89e4eSSatish Balay 
588e5c89e4eSSatish Balay    Input Parameters:
589e5c89e4eSSatish Balay +  n  - number of values
59059e16d97SJunchao Zhang .  X  - array of integers
59159e16d97SJunchao Zhang -  Y  - second array of scalars
592e5c89e4eSSatish Balay 
593e5c89e4eSSatish Balay    Level: intermediate
594e5c89e4eSSatish Balay 
595e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
596e5c89e4eSSatish Balay @*/
597b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithScalarArray(PetscInt n,PetscInt X[],PetscScalar Y[])
598e5c89e4eSSatish Balay {
599e5c89e4eSSatish Balay   PetscErrorCode ierr;
60059e16d97SJunchao Zhang   PetscInt       pivot,t1;
60159e16d97SJunchao Zhang   PetscScalar    t2;
602e5c89e4eSSatish Balay 
603e5c89e4eSSatish Balay   PetscFunctionBegin;
60459e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithScalarArray,X,Y,n,pivot,t1,t2,ierr);
60517df18f2SToby Isaac   PetscFunctionReturn(0);
60617df18f2SToby Isaac }
60717df18f2SToby Isaac 
6086c2863d0SBarry Smith /*@C
60917df18f2SToby Isaac    PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order;
61017df18f2SToby Isaac        changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
61117df18f2SToby Isaac        provide workspace (the size of an element in the data array) to use when sorting.
61217df18f2SToby Isaac 
61317df18f2SToby Isaac    Not Collective
61417df18f2SToby Isaac 
61517df18f2SToby Isaac    Input Parameters:
61617df18f2SToby Isaac +  n  - number of values
61759e16d97SJunchao Zhang .  X  - array of integers
61859e16d97SJunchao Zhang .  Y  - second array of data
61917df18f2SToby Isaac .  size - sizeof elements in the data array in bytes
62059e16d97SJunchao Zhang -  t2   - workspace of "size" bytes used when sorting
62117df18f2SToby Isaac 
62217df18f2SToby Isaac    Level: intermediate
62317df18f2SToby Isaac 
62417df18f2SToby Isaac .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
62517df18f2SToby Isaac @*/
626b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithDataArray(PetscInt n,PetscInt X[],void *Y,size_t size,void *t2)
62717df18f2SToby Isaac {
62817df18f2SToby Isaac   PetscErrorCode ierr;
62959e16d97SJunchao Zhang   char           *YY = (char*)Y;
63059e16d97SJunchao Zhang   PetscInt       i,j,p,t1,pivot,hi=n-1,l,r;
63117df18f2SToby Isaac 
63217df18f2SToby Isaac   PetscFunctionBegin;
63317df18f2SToby Isaac   if (n<8) {
63459e16d97SJunchao Zhang     for (i=0; i<n; i++) {
63559e16d97SJunchao Zhang       pivot = X[i];
63659e16d97SJunchao Zhang       for (j=i+1; j<n; j++) {
63759e16d97SJunchao Zhang         if (pivot > X[j]) {
63859e16d97SJunchao Zhang           SWAP2Data(X[i],X[j],YY+size*i,YY+size*j,t1,t2,size);
63959e16d97SJunchao Zhang           pivot = X[i];
64017df18f2SToby Isaac         }
64117df18f2SToby Isaac       }
64217df18f2SToby Isaac     }
64317df18f2SToby Isaac   } else {
64459e16d97SJunchao Zhang     /* Two way partition */
64559e16d97SJunchao Zhang     p     = MEDIAN(X,hi);
64659e16d97SJunchao Zhang     pivot = X[p];
64759e16d97SJunchao Zhang     l     = 0;
64859e16d97SJunchao Zhang     r     = hi;
64959e16d97SJunchao Zhang     while(1) {
65059e16d97SJunchao Zhang       while (X[l] < pivot) l++;
65159e16d97SJunchao Zhang       while (X[r] > pivot) r--;
65259e16d97SJunchao Zhang       if (l >= r) {r++; break;}
65359e16d97SJunchao Zhang       SWAP2Data(X[l],X[r],YY+size*l,YY+size*r,t1,t2,size);
65459e16d97SJunchao Zhang       l++;
65559e16d97SJunchao Zhang       r--;
65659e16d97SJunchao Zhang     }
65759e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(l,X,Y,size,t2);CHKERRQ(ierr);
65859e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(hi-r+1,X+r,YY+size*r,size,t2);CHKERRQ(ierr);
65917df18f2SToby Isaac   }
66017df18f2SToby Isaac   PetscFunctionReturn(0);
66117df18f2SToby Isaac }
66217df18f2SToby Isaac 
66321e72a00SBarry Smith /*@
66421e72a00SBarry Smith    PetscMergeIntArray -     Merges two SORTED integer arrays, removes duplicate elements.
66521e72a00SBarry Smith 
66621e72a00SBarry Smith    Not Collective
66721e72a00SBarry Smith 
66821e72a00SBarry Smith    Input Parameters:
66921e72a00SBarry Smith +  an  - number of values in the first array
67021e72a00SBarry Smith .  aI  - first sorted array of integers
67121e72a00SBarry Smith .  bn  - number of values in the second array
67221e72a00SBarry Smith -  bI  - second array of integers
67321e72a00SBarry Smith 
67421e72a00SBarry Smith    Output Parameters:
67521e72a00SBarry Smith +  n   - number of values in the merged array
6766c2863d0SBarry Smith -  L   - merged sorted array, this is allocated if an array is not provided
67721e72a00SBarry Smith 
67821e72a00SBarry Smith    Level: intermediate
67921e72a00SBarry Smith 
68021e72a00SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
68121e72a00SBarry Smith @*/
6826c2863d0SBarry Smith PetscErrorCode  PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
68321e72a00SBarry Smith {
68421e72a00SBarry Smith   PetscErrorCode ierr;
68521e72a00SBarry Smith   PetscInt       *L_ = *L, ak, bk, k;
68621e72a00SBarry Smith 
68721e72a00SBarry Smith   if (!L_) {
68821e72a00SBarry Smith     ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr);
68921e72a00SBarry Smith     L_   = *L;
69021e72a00SBarry Smith   }
69121e72a00SBarry Smith   k = ak = bk = 0;
69221e72a00SBarry Smith   while (ak < an && bk < bn) {
69321e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
69421e72a00SBarry Smith       L_[k] = aI[ak];
69521e72a00SBarry Smith       ++ak;
69621e72a00SBarry Smith       ++bk;
69721e72a00SBarry Smith       ++k;
69821e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
69921e72a00SBarry Smith       L_[k] = aI[ak];
70021e72a00SBarry Smith       ++ak;
70121e72a00SBarry Smith       ++k;
70221e72a00SBarry Smith     } else {
70321e72a00SBarry Smith       L_[k] = bI[bk];
70421e72a00SBarry Smith       ++bk;
70521e72a00SBarry Smith       ++k;
70621e72a00SBarry Smith     }
70721e72a00SBarry Smith   }
70821e72a00SBarry Smith   if (ak < an) {
709580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
71021e72a00SBarry Smith     k   += (an-ak);
71121e72a00SBarry Smith   }
71221e72a00SBarry Smith   if (bk < bn) {
713580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
71421e72a00SBarry Smith     k   += (bn-bk);
71521e72a00SBarry Smith   }
71621e72a00SBarry Smith   *n = k;
71721e72a00SBarry Smith   PetscFunctionReturn(0);
71821e72a00SBarry Smith }
719b4301105SBarry Smith 
720c1f0200aSDmitry Karpeev /*@
72121e72a00SBarry Smith    PetscMergeIntArrayPair -     Merges two SORTED integer arrays that share NO common values along with an additional array of integers.
722c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
723c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
724c1f0200aSDmitry Karpeev 
725c1f0200aSDmitry Karpeev    Not Collective
726c1f0200aSDmitry Karpeev 
727c1f0200aSDmitry Karpeev    Input Parameters:
728c1f0200aSDmitry Karpeev +  an  - number of values in the first array
729c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
730c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
731c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
732c1f0200aSDmitry Karpeev .  bI  - second array of integers
733c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
734c1f0200aSDmitry Karpeev 
735c1f0200aSDmitry Karpeev    Output Parameters:
736c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
73714c49006SJed Brown .  L   - merged sorted array
738c1f0200aSDmitry Karpeev -  J   - merged additional array
739c1f0200aSDmitry Karpeev 
74095452b02SPatrick Sanan    Notes:
74195452b02SPatrick 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
742c1f0200aSDmitry Karpeev    Level: intermediate
743c1f0200aSDmitry Karpeev 
744c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
745c1f0200aSDmitry Karpeev @*/
7466c2863d0SBarry 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)
747c1f0200aSDmitry Karpeev {
748c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
74970d8d27cSBarry Smith   PetscInt       n_, *L_, *J_, ak, bk, k;
750c1f0200aSDmitry Karpeev 
75114c49006SJed Brown   PetscFunctionBegin;
75270d8d27cSBarry Smith   PetscValidIntPointer(L,8);
75370d8d27cSBarry Smith   PetscValidIntPointer(J,9);
754c1f0200aSDmitry Karpeev   n_ = an + bn;
755c1f0200aSDmitry Karpeev   *n = n_;
75670d8d27cSBarry Smith   if (!*L) {
757785e854fSJed Brown     ierr = PetscMalloc1(n_, L);CHKERRQ(ierr);
75870d8d27cSBarry Smith   }
759d7aa01a8SHong Zhang   L_ = *L;
76070d8d27cSBarry Smith   if (!*J) {
76170d8d27cSBarry Smith     ierr = PetscMalloc1(n_, J);CHKERRQ(ierr);
762c1f0200aSDmitry Karpeev   }
763c1f0200aSDmitry Karpeev   J_   = *J;
764c1f0200aSDmitry Karpeev   k = ak = bk = 0;
765c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
766c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
767d7aa01a8SHong Zhang       L_[k] = aI[ak];
768c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
769c1f0200aSDmitry Karpeev       ++ak;
770c1f0200aSDmitry Karpeev       ++k;
771a297a907SKarl Rupp     } else {
772d7aa01a8SHong Zhang       L_[k] = bI[bk];
773c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
774c1f0200aSDmitry Karpeev       ++bk;
775c1f0200aSDmitry Karpeev       ++k;
776c1f0200aSDmitry Karpeev     }
777c1f0200aSDmitry Karpeev   }
778c1f0200aSDmitry Karpeev   if (ak < an) {
779580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
780580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,aJ+ak,an-ak);CHKERRQ(ierr);
781c1f0200aSDmitry Karpeev     k   += (an-ak);
782c1f0200aSDmitry Karpeev   }
783c1f0200aSDmitry Karpeev   if (bk < bn) {
784580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
785580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,bJ+bk,bn-bk);CHKERRQ(ierr);
786c1f0200aSDmitry Karpeev   }
787c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
788c1f0200aSDmitry Karpeev }
789c1f0200aSDmitry Karpeev 
790e498c390SJed Brown /*@
791e498c390SJed Brown    PetscMergeMPIIntArray -     Merges two SORTED integer arrays.
792e498c390SJed Brown 
793e498c390SJed Brown    Not Collective
794e498c390SJed Brown 
795e498c390SJed Brown    Input Parameters:
796e498c390SJed Brown +  an  - number of values in the first array
797e498c390SJed Brown .  aI  - first sorted array of integers
798e498c390SJed Brown .  bn  - number of values in the second array
799e498c390SJed Brown -  bI  - second array of integers
800e498c390SJed Brown 
801e498c390SJed Brown    Output Parameters:
802e498c390SJed Brown +  n   - number of values in the merged array (<= an + bn)
803e498c390SJed Brown -  L   - merged sorted array, allocated if address of NULL pointer is passed
804e498c390SJed Brown 
805e498c390SJed Brown    Level: intermediate
806e498c390SJed Brown 
807e498c390SJed Brown .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
808e498c390SJed Brown @*/
809e498c390SJed Brown PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L)
810e498c390SJed Brown {
811e498c390SJed Brown   PetscErrorCode ierr;
812e498c390SJed Brown   PetscInt       ai,bi,k;
813e498c390SJed Brown 
814e498c390SJed Brown   PetscFunctionBegin;
815e498c390SJed Brown   if (!*L) {ierr = PetscMalloc1((an+bn),L);CHKERRQ(ierr);}
816e498c390SJed Brown   for (ai=0,bi=0,k=0; ai<an || bi<bn; ) {
817e498c390SJed Brown     PetscInt t = -1;
818e498c390SJed Brown     for ( ; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
819e498c390SJed Brown     for ( ; bi<bn && bI[bi] == t; bi++);
820e498c390SJed Brown     for ( ; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
821e498c390SJed Brown     for ( ; ai<an && aI[ai] == t; ai++);
822e498c390SJed Brown   }
823e498c390SJed Brown   *n = k;
824e498c390SJed Brown   PetscFunctionReturn(0);
825e498c390SJed Brown }
826e5c89e4eSSatish Balay 
8276c2863d0SBarry Smith /*@C
82842eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
82942eaa7f3SBarry Smith 
83042eaa7f3SBarry Smith    Not Collective
83142eaa7f3SBarry Smith 
83242eaa7f3SBarry Smith    Input Parameters:
83342eaa7f3SBarry Smith +  n  - number of values
83442eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
83542eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
83642eaa7f3SBarry Smith 
83742eaa7f3SBarry Smith    Output Parameters:
83842eaa7f3SBarry Smith +  Nlevels - the number of levels
83942eaa7f3SBarry Smith .  Level - for each node tells its level
84042eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
84142eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
84242eaa7f3SBarry Smith -  Column - for each id tells its column index
84342eaa7f3SBarry Smith 
8446c2863d0SBarry Smith    Level: developer
84542eaa7f3SBarry Smith 
84695452b02SPatrick Sanan    Notes:
84795452b02SPatrick Sanan     This code is not currently used
84842eaa7f3SBarry Smith 
84942eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation()
85042eaa7f3SBarry Smith @*/
8517087cfbeSBarry Smith PetscErrorCode  PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column)
85242eaa7f3SBarry Smith {
85342eaa7f3SBarry Smith   PetscInt       i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column;
85442eaa7f3SBarry Smith   PetscErrorCode ierr;
855ace3abfcSBarry Smith   PetscBool      done = PETSC_FALSE;
85642eaa7f3SBarry Smith 
85742eaa7f3SBarry Smith   PetscFunctionBegin;
85842eaa7f3SBarry Smith   if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set");
85942eaa7f3SBarry Smith   for (i=0; i<n; i++) {
86042eaa7f3SBarry Smith     if (mask[i]) continue;
86142eaa7f3SBarry Smith     if (parentid[i]  == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent");
86242eaa7f3SBarry Smith     if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked");
86342eaa7f3SBarry Smith   }
86442eaa7f3SBarry Smith 
86542eaa7f3SBarry Smith   for (i=0; i<n; i++) {
86642eaa7f3SBarry Smith     if (!mask[i]) nmask++;
86742eaa7f3SBarry Smith   }
86842eaa7f3SBarry Smith 
86942eaa7f3SBarry Smith   /* determine the level in the tree of each node */
8701795a4d1SJed Brown   ierr = PetscCalloc1(n,&level);CHKERRQ(ierr);
871a297a907SKarl Rupp 
87242eaa7f3SBarry Smith   level[0] = 1;
87342eaa7f3SBarry Smith   while (!done) {
87442eaa7f3SBarry Smith     done = PETSC_TRUE;
87542eaa7f3SBarry Smith     for (i=0; i<n; i++) {
87642eaa7f3SBarry Smith       if (mask[i]) continue;
87742eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
87842eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
87942eaa7f3SBarry Smith     }
88042eaa7f3SBarry Smith   }
88142eaa7f3SBarry Smith   for (i=0; i<n; i++) {
88242eaa7f3SBarry Smith     level[i]--;
88342eaa7f3SBarry Smith     nlevels = PetscMax(nlevels,level[i]);
88442eaa7f3SBarry Smith   }
88542eaa7f3SBarry Smith 
88642eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
8871795a4d1SJed Brown   ierr = PetscCalloc1(nlevels,&levelcnt);CHKERRQ(ierr);
88842eaa7f3SBarry Smith   for (i=0; i<n; i++) {
88942eaa7f3SBarry Smith     if (mask[i]) continue;
89042eaa7f3SBarry Smith     levelcnt[level[i]-1]++;
89142eaa7f3SBarry Smith   }
892a297a907SKarl Rupp   for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]);
89342eaa7f3SBarry Smith 
89442eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
895dcca6d9dSJed Brown   ierr = PetscMalloc2(levelmax,&workid,levelmax,&workparentid);CHKERRQ(ierr);
896785e854fSJed Brown   ierr = PetscMalloc1(nmask,&idbylevel);CHKERRQ(ierr);
89742eaa7f3SBarry Smith   for (j=1; j<=nlevels;j++) {
89842eaa7f3SBarry Smith     cnt = 0;
89942eaa7f3SBarry Smith     for (i=0; i<n; i++) {
90042eaa7f3SBarry Smith       if (mask[i]) continue;
90142eaa7f3SBarry Smith       if (level[i] != j) continue;
90242eaa7f3SBarry Smith       workid[cnt]         = i;
90342eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
90442eaa7f3SBarry Smith     }
90542eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
90642eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
90742eaa7f3SBarry Smith     ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr);
90842eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
90942eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
910580bdb30SBarry Smith     ierr  = PetscArraycpy(idbylevel+tcnt,workid,cnt);CHKERRQ(ierr);
91142eaa7f3SBarry Smith     tcnt += cnt;
91242eaa7f3SBarry Smith   }
91342eaa7f3SBarry Smith   if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes");
91442eaa7f3SBarry Smith   ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr);
91542eaa7f3SBarry Smith 
91642eaa7f3SBarry Smith   /* for each node list its column */
917785e854fSJed Brown   ierr = PetscMalloc1(n,&column);CHKERRQ(ierr);
91842eaa7f3SBarry Smith   cnt = 0;
91942eaa7f3SBarry Smith   for (j=0; j<nlevels; j++) {
92042eaa7f3SBarry Smith     for (i=0; i<levelcnt[j]; i++) {
92142eaa7f3SBarry Smith       column[idbylevel[cnt++]] = i;
92242eaa7f3SBarry Smith     }
92342eaa7f3SBarry Smith   }
92442eaa7f3SBarry Smith 
92542eaa7f3SBarry Smith   *Nlevels   = nlevels;
92642eaa7f3SBarry Smith   *Level     = level;
92742eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
92842eaa7f3SBarry Smith   *Idbylevel = idbylevel;
92942eaa7f3SBarry Smith   *Column    = column;
93042eaa7f3SBarry Smith   PetscFunctionReturn(0);
93142eaa7f3SBarry Smith }
932ce605777SToby Isaac 
933ce605777SToby Isaac /*@
934ce605777SToby Isaac   PetscParallelSortedInt - Check whether an integer array, distributed over a communicator, is globally sorted.
935ce605777SToby Isaac 
936ce605777SToby Isaac   Collective
937ce605777SToby Isaac 
938ce605777SToby Isaac   Input Parameters:
939ce605777SToby Isaac + comm - the MPI communicator
940ce605777SToby Isaac . n - the local number of integers
941ce605777SToby Isaac - keys - the local array of integers
942ce605777SToby Isaac 
943ce605777SToby Isaac   Output Parameters:
944ce605777SToby Isaac . is_sorted - whether the array is globally sorted
945ce605777SToby Isaac 
946ce605777SToby Isaac   Level: developer
947ce605777SToby Isaac 
948ce605777SToby Isaac .seealso: PetscParallelSortInt()
949ce605777SToby Isaac @*/
950ce605777SToby Isaac PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
951ce605777SToby Isaac {
952ce605777SToby Isaac   PetscBool      sorted;
953ce605777SToby Isaac   PetscInt       i, min, max, prevmax;
954*2a1da528SToby Isaac   PetscMPIInt    rank;
955ce605777SToby Isaac   PetscErrorCode ierr;
956ce605777SToby Isaac 
957ce605777SToby Isaac   PetscFunctionBegin;
958ce605777SToby Isaac   sorted = PETSC_TRUE;
959ce605777SToby Isaac   min    = PETSC_MAX_INT;
960ce605777SToby Isaac   max    = PETSC_MIN_INT;
961ce605777SToby Isaac   if (n) {
962ce605777SToby Isaac     min = keys[0];
963ce605777SToby Isaac     max = keys[0];
964ce605777SToby Isaac   }
965ce605777SToby Isaac   for (i = 1; i < n; i++) {
966ce605777SToby Isaac     if (keys[i] < keys[i - 1]) break;
967ce605777SToby Isaac     min = PetscMin(min,keys[i]);
968ce605777SToby Isaac     max = PetscMax(max,keys[i]);
969ce605777SToby Isaac   }
970ce605777SToby Isaac   if (i < n) sorted = PETSC_FALSE;
971ce605777SToby Isaac   prevmax = PETSC_MIN_INT;
972ce605777SToby Isaac   ierr = MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
973*2a1da528SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
974*2a1da528SToby Isaac   if (!rank) prevmax = PETSC_MIN_INT;
975ce605777SToby Isaac   if (prevmax > min) sorted = PETSC_FALSE;
976ce605777SToby Isaac   ierr = MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr);
977ce605777SToby Isaac   PetscFunctionReturn(0);
978ce605777SToby Isaac }
979ce605777SToby Isaac 
980