xref: /petsc/src/sys/utils/sorti.c (revision fb61b9e4c10773d1e08e7a826aed0c51ec8ec9ed)
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 /*@
2146a7c706bSVaclav Hapla    PetscSortedInt - Determines whether the array is sorted.
2156a7c706bSVaclav Hapla 
2166a7c706bSVaclav Hapla    Not Collective
2176a7c706bSVaclav Hapla 
2186a7c706bSVaclav Hapla    Input Parameters:
2196a7c706bSVaclav Hapla +  n  - number of values
2206a7c706bSVaclav Hapla -  X  - array of integers
2216a7c706bSVaclav Hapla 
2226a7c706bSVaclav Hapla    Output Parameters:
2236a7c706bSVaclav Hapla .  sorted - flag whether the array is sorted
2246a7c706bSVaclav Hapla 
2256a7c706bSVaclav Hapla    Level: intermediate
2266a7c706bSVaclav Hapla 
2276a7c706bSVaclav Hapla .seealso: PetscSortInt(), PetscSortedMPIInt(), PetscSortedReal()
2286a7c706bSVaclav Hapla @*/
2296a7c706bSVaclav Hapla PetscErrorCode  PetscSortedInt(PetscInt n,const PetscInt X[],PetscBool *sorted)
2306a7c706bSVaclav Hapla {
2316a7c706bSVaclav Hapla   PetscFunctionBegin;
2326a7c706bSVaclav Hapla   PetscSorted(n,X,*sorted);
2336a7c706bSVaclav Hapla   PetscFunctionReturn(0);
2346a7c706bSVaclav Hapla }
2356a7c706bSVaclav Hapla 
2366a7c706bSVaclav Hapla /*@
237e5c89e4eSSatish Balay    PetscSortInt - Sorts an array of integers in place in increasing order.
238e5c89e4eSSatish Balay 
239e5c89e4eSSatish Balay    Not Collective
240e5c89e4eSSatish Balay 
241e5c89e4eSSatish Balay    Input Parameters:
242e5c89e4eSSatish Balay +  n  - number of values
24359e16d97SJunchao Zhang -  X  - array of integers
244e5c89e4eSSatish Balay 
245e5c89e4eSSatish Balay    Level: intermediate
246e5c89e4eSSatish Balay 
247e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntWithPermutation()
248e5c89e4eSSatish Balay @*/
249b4f531b9SJunchao Zhang PetscErrorCode  PetscSortInt(PetscInt n,PetscInt X[])
250e5c89e4eSSatish Balay {
25159e16d97SJunchao Zhang   PetscErrorCode ierr;
25259e16d97SJunchao Zhang   PetscInt       pivot,t1;
253e5c89e4eSSatish Balay 
254e5c89e4eSSatish Balay   PetscFunctionBegin;
25559e16d97SJunchao Zhang   QuickSort1(PetscSortInt,X,n,pivot,t1,ierr);
256e5c89e4eSSatish Balay   PetscFunctionReturn(0);
257e5c89e4eSSatish Balay }
258e5c89e4eSSatish Balay 
2591db36a52SBarry Smith /*@
260ce605777SToby Isaac    PetscSortReverseInt - Sorts an array of integers in place in decreasing order.
261ce605777SToby Isaac 
262ce605777SToby Isaac    Not Collective
263ce605777SToby Isaac 
264ce605777SToby Isaac    Input Parameters:
265ce605777SToby Isaac +  n  - number of values
266ce605777SToby Isaac -  X  - array of integers
267ce605777SToby Isaac 
268ce605777SToby Isaac    Level: intermediate
269ce605777SToby Isaac 
270ce605777SToby Isaac .seealso: PetscSortInt(), PetscSortIntWithPermutation()
271ce605777SToby Isaac @*/
272ce605777SToby Isaac PetscErrorCode  PetscSortReverseInt(PetscInt n,PetscInt X[])
273ce605777SToby Isaac {
274ce605777SToby Isaac   PetscErrorCode ierr;
275ce605777SToby Isaac   PetscInt       pivot,t1;
276ce605777SToby Isaac 
277ce605777SToby Isaac   PetscFunctionBegin;
278ce605777SToby Isaac   QuickSortReverse1(PetscSortReverseInt,X,n,pivot,t1,ierr);
279ce605777SToby Isaac   PetscFunctionReturn(0);
280ce605777SToby Isaac }
281ce605777SToby Isaac 
282ce605777SToby Isaac /*@
28322ab5688SLisandro Dalcin    PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted input array
28422ab5688SLisandro Dalcin 
28522ab5688SLisandro Dalcin    Not Collective
28622ab5688SLisandro Dalcin 
28722ab5688SLisandro Dalcin    Input Parameters:
28822ab5688SLisandro Dalcin +  n  - number of values
28959e16d97SJunchao Zhang -  X  - sorted array of integers
29022ab5688SLisandro Dalcin 
29122ab5688SLisandro Dalcin    Output Parameter:
29222ab5688SLisandro Dalcin .  n - number of non-redundant values
29322ab5688SLisandro Dalcin 
29422ab5688SLisandro Dalcin    Level: intermediate
29522ab5688SLisandro Dalcin 
29622ab5688SLisandro Dalcin .seealso: PetscSortInt()
29722ab5688SLisandro Dalcin @*/
29859e16d97SJunchao Zhang PetscErrorCode  PetscSortedRemoveDupsInt(PetscInt *n,PetscInt X[])
29922ab5688SLisandro Dalcin {
30022ab5688SLisandro Dalcin   PetscInt i,s = 0,N = *n, b = 0;
30122ab5688SLisandro Dalcin 
30222ab5688SLisandro Dalcin   PetscFunctionBegin;
303*fb61b9e4SStefano Zampini   PetscCheckSorted(*n,X);
30422ab5688SLisandro Dalcin   for (i=0; i<N-1; i++) {
30559e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
30659e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
30722ab5688SLisandro Dalcin     } else s++;
30822ab5688SLisandro Dalcin   }
30922ab5688SLisandro Dalcin   *n = N - s;
31022ab5688SLisandro Dalcin   PetscFunctionReturn(0);
31122ab5688SLisandro Dalcin }
31222ab5688SLisandro Dalcin 
31322ab5688SLisandro Dalcin /*@
3141db36a52SBarry Smith    PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries
3151db36a52SBarry Smith 
3161db36a52SBarry Smith    Not Collective
3171db36a52SBarry Smith 
3181db36a52SBarry Smith    Input Parameters:
3191db36a52SBarry Smith +  n  - number of values
32059e16d97SJunchao Zhang -  X  - array of integers
3211db36a52SBarry Smith 
3221db36a52SBarry Smith    Output Parameter:
3231db36a52SBarry Smith .  n - number of non-redundant values
3241db36a52SBarry Smith 
3251db36a52SBarry Smith    Level: intermediate
3261db36a52SBarry Smith 
32722ab5688SLisandro Dalcin .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortedRemoveDupsInt()
3281db36a52SBarry Smith @*/
32959e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsInt(PetscInt *n,PetscInt X[])
3301db36a52SBarry Smith {
3311db36a52SBarry Smith   PetscErrorCode ierr;
3321db36a52SBarry Smith 
3331db36a52SBarry Smith   PetscFunctionBegin;
33459e16d97SJunchao Zhang   ierr = PetscSortInt(*n,X);CHKERRQ(ierr);
33559e16d97SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(n,X);CHKERRQ(ierr);
3361db36a52SBarry Smith   PetscFunctionReturn(0);
3371db36a52SBarry Smith }
3381db36a52SBarry Smith 
33960e03357SMatthew G Knepley /*@
340d9f1162dSJed Brown   PetscFindInt - Finds integer in a sorted array of integers
34160e03357SMatthew G Knepley 
34260e03357SMatthew G Knepley    Not Collective
34360e03357SMatthew G Knepley 
34460e03357SMatthew G Knepley    Input Parameters:
34560e03357SMatthew G Knepley +  key - the integer to locate
346d9f1162dSJed Brown .  n   - number of values in the array
34759e16d97SJunchao Zhang -  X  - array of integers
34860e03357SMatthew G Knepley 
34960e03357SMatthew G Knepley    Output Parameter:
350d9f1162dSJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
35160e03357SMatthew G Knepley 
35260e03357SMatthew G Knepley    Level: intermediate
35360e03357SMatthew G Knepley 
35460e03357SMatthew G Knepley .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
35560e03357SMatthew G Knepley @*/
35659e16d97SJunchao Zhang PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
35760e03357SMatthew G Knepley {
358d9f1162dSJed Brown   PetscInt lo = 0,hi = n;
35960e03357SMatthew G Knepley 
36060e03357SMatthew G Knepley   PetscFunctionBegin;
36160e03357SMatthew G Knepley   PetscValidPointer(loc,4);
36298781d82SMatthew G Knepley   if (!n) {*loc = -1; PetscFunctionReturn(0);}
36359e16d97SJunchao Zhang   PetscValidPointer(X,3);
3646a7c706bSVaclav Hapla   PetscCheckSorted(n,X);
365d9f1162dSJed Brown   while (hi - lo > 1) {
366d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo)/2;
36759e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
368d9f1162dSJed Brown     else               lo = mid;
36960e03357SMatthew G Knepley   }
37059e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
37160e03357SMatthew G Knepley   PetscFunctionReturn(0);
37260e03357SMatthew G Knepley }
37360e03357SMatthew G Knepley 
374d2aeb606SJed Brown /*@
375f1cab4e1SJunchao Zhang   PetscCheckDupsInt - Checks if an integer array has duplicates
376f1cab4e1SJunchao Zhang 
377f1cab4e1SJunchao Zhang    Not Collective
378f1cab4e1SJunchao Zhang 
379f1cab4e1SJunchao Zhang    Input Parameters:
380f1cab4e1SJunchao Zhang +  n  - number of values in the array
381f1cab4e1SJunchao Zhang -  X  - array of integers
382f1cab4e1SJunchao Zhang 
383f1cab4e1SJunchao Zhang 
384f1cab4e1SJunchao Zhang    Output Parameter:
385f1cab4e1SJunchao Zhang .  dups - True if the array has dups, otherwise false
386f1cab4e1SJunchao Zhang 
387f1cab4e1SJunchao Zhang    Level: intermediate
388f1cab4e1SJunchao Zhang 
389f1cab4e1SJunchao Zhang .seealso: PetscSortRemoveDupsInt()
390f1cab4e1SJunchao Zhang @*/
391f1cab4e1SJunchao Zhang PetscErrorCode PetscCheckDupsInt(PetscInt n,const PetscInt X[],PetscBool *dups)
392f1cab4e1SJunchao Zhang {
393f1cab4e1SJunchao Zhang   PetscErrorCode ierr;
394f1cab4e1SJunchao Zhang   PetscInt       i;
395f1cab4e1SJunchao Zhang   PetscHSetI     ht;
396f1cab4e1SJunchao Zhang   PetscBool      missing;
397f1cab4e1SJunchao Zhang 
398f1cab4e1SJunchao Zhang   PetscFunctionBegin;
399f1cab4e1SJunchao Zhang   PetscValidPointer(dups,3);
400f1cab4e1SJunchao Zhang   *dups = PETSC_FALSE;
401f1cab4e1SJunchao Zhang   if (n > 1) {
402f1cab4e1SJunchao Zhang     ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
403f1cab4e1SJunchao Zhang     ierr = PetscHSetIResize(ht,n);CHKERRQ(ierr);
404f1cab4e1SJunchao Zhang     for (i=0; i<n; i++) {
405f1cab4e1SJunchao Zhang       ierr = PetscHSetIQueryAdd(ht,X[i],&missing);CHKERRQ(ierr);
406f1cab4e1SJunchao Zhang       if(!missing) {*dups = PETSC_TRUE; break;}
407f1cab4e1SJunchao Zhang     }
408f1cab4e1SJunchao Zhang     ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
409f1cab4e1SJunchao Zhang   }
410f1cab4e1SJunchao Zhang   PetscFunctionReturn(0);
411f1cab4e1SJunchao Zhang }
412f1cab4e1SJunchao Zhang 
413f1cab4e1SJunchao Zhang /*@
414d2aeb606SJed Brown   PetscFindMPIInt - Finds MPI integer in a sorted array of integers
415d2aeb606SJed Brown 
416d2aeb606SJed Brown    Not Collective
417d2aeb606SJed Brown 
418d2aeb606SJed Brown    Input Parameters:
419d2aeb606SJed Brown +  key - the integer to locate
420d2aeb606SJed Brown .  n   - number of values in the array
42159e16d97SJunchao Zhang -  X   - array of integers
422d2aeb606SJed Brown 
423d2aeb606SJed Brown    Output Parameter:
424d2aeb606SJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
425d2aeb606SJed Brown 
426d2aeb606SJed Brown    Level: intermediate
427d2aeb606SJed Brown 
428d2aeb606SJed Brown .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
429d2aeb606SJed Brown @*/
43059e16d97SJunchao Zhang PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
431d2aeb606SJed Brown {
432d2aeb606SJed Brown   PetscInt lo = 0,hi = n;
433d2aeb606SJed Brown 
434d2aeb606SJed Brown   PetscFunctionBegin;
435d2aeb606SJed Brown   PetscValidPointer(loc,4);
436d2aeb606SJed Brown   if (!n) {*loc = -1; PetscFunctionReturn(0);}
43759e16d97SJunchao Zhang   PetscValidPointer(X,3);
4386a7c706bSVaclav Hapla   PetscCheckSorted(n,X);
439d2aeb606SJed Brown   while (hi - lo > 1) {
440d2aeb606SJed Brown     PetscInt mid = lo + (hi - lo)/2;
44159e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
442d2aeb606SJed Brown     else               lo = mid;
443d2aeb606SJed Brown   }
44459e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
445e5c89e4eSSatish Balay   PetscFunctionReturn(0);
446e5c89e4eSSatish Balay }
447e5c89e4eSSatish Balay 
448e5c89e4eSSatish Balay /*@
449e5c89e4eSSatish Balay    PetscSortIntWithArray - Sorts an array of integers in place in increasing order;
450e5c89e4eSSatish Balay        changes a second array to match the sorted first array.
451e5c89e4eSSatish Balay 
452e5c89e4eSSatish Balay    Not Collective
453e5c89e4eSSatish Balay 
454e5c89e4eSSatish Balay    Input Parameters:
455e5c89e4eSSatish Balay +  n  - number of values
45659e16d97SJunchao Zhang .  X  - array of integers
45759e16d97SJunchao Zhang -  Y  - second array of integers
458e5c89e4eSSatish Balay 
459e5c89e4eSSatish Balay    Level: intermediate
460e5c89e4eSSatish Balay 
461e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
462e5c89e4eSSatish Balay @*/
463b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithArray(PetscInt n,PetscInt X[],PetscInt Y[])
464e5c89e4eSSatish Balay {
465e5c89e4eSSatish Balay   PetscErrorCode ierr;
46659e16d97SJunchao Zhang   PetscInt       pivot,t1,t2;
467e5c89e4eSSatish Balay 
468e5c89e4eSSatish Balay   PetscFunctionBegin;
46959e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithArray,X,Y,n,pivot,t1,t2,ierr);
470c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
471c1f0200aSDmitry Karpeev }
472c1f0200aSDmitry Karpeev 
473c1f0200aSDmitry Karpeev /*@
474c1f0200aSDmitry Karpeev    PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order;
475c1f0200aSDmitry Karpeev        changes a pair of integer arrays to match the sorted first array.
476c1f0200aSDmitry Karpeev 
477c1f0200aSDmitry Karpeev    Not Collective
478c1f0200aSDmitry Karpeev 
479c1f0200aSDmitry Karpeev    Input Parameters:
480c1f0200aSDmitry Karpeev +  n  - number of values
48159e16d97SJunchao Zhang .  X  - array of integers
48259e16d97SJunchao Zhang .  Y  - second array of integers (first array of the pair)
48359e16d97SJunchao Zhang -  Z  - third array of integers  (second array of the pair)
484c1f0200aSDmitry Karpeev 
485c1f0200aSDmitry Karpeev    Level: intermediate
486c1f0200aSDmitry Karpeev 
487c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray()
488c1f0200aSDmitry Karpeev @*/
489b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithArrayPair(PetscInt n,PetscInt X[],PetscInt Y[],PetscInt Z[])
490c1f0200aSDmitry Karpeev {
491c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
49259e16d97SJunchao Zhang   PetscInt       pivot,t1,t2,t3;
493c1f0200aSDmitry Karpeev 
494c1f0200aSDmitry Karpeev   PetscFunctionBegin;
49559e16d97SJunchao Zhang   QuickSort3(PetscSortIntWithArrayPair,X,Y,Z,n,pivot,t1,t2,t3,ierr);
496c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
497c1f0200aSDmitry Karpeev }
498c1f0200aSDmitry Karpeev 
49917d7d925SStefano Zampini /*@
5006a7c706bSVaclav Hapla    PetscSortedMPIInt - Determines whether the array is sorted.
5016a7c706bSVaclav Hapla 
5026a7c706bSVaclav Hapla    Not Collective
5036a7c706bSVaclav Hapla 
5046a7c706bSVaclav Hapla    Input Parameters:
5056a7c706bSVaclav Hapla +  n  - number of values
5066a7c706bSVaclav Hapla -  X  - array of integers
5076a7c706bSVaclav Hapla 
5086a7c706bSVaclav Hapla    Output Parameters:
5096a7c706bSVaclav Hapla .  sorted - flag whether the array is sorted
5106a7c706bSVaclav Hapla 
5116a7c706bSVaclav Hapla    Level: intermediate
5126a7c706bSVaclav Hapla 
5136a7c706bSVaclav Hapla .seealso: PetscSortMPIInt(), PetscSortedInt(), PetscSortedReal()
5146a7c706bSVaclav Hapla @*/
5156a7c706bSVaclav Hapla PetscErrorCode  PetscSortedMPIInt(PetscInt n,const PetscMPIInt X[],PetscBool *sorted)
5166a7c706bSVaclav Hapla {
5176a7c706bSVaclav Hapla   PetscFunctionBegin;
5186a7c706bSVaclav Hapla   PetscSorted(n,X,*sorted);
5196a7c706bSVaclav Hapla   PetscFunctionReturn(0);
5206a7c706bSVaclav Hapla }
5216a7c706bSVaclav Hapla 
5226a7c706bSVaclav Hapla /*@
52317d7d925SStefano Zampini    PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order.
52417d7d925SStefano Zampini 
52517d7d925SStefano Zampini    Not Collective
52617d7d925SStefano Zampini 
52717d7d925SStefano Zampini    Input Parameters:
52817d7d925SStefano Zampini +  n  - number of values
52959e16d97SJunchao Zhang -  X  - array of integers
53017d7d925SStefano Zampini 
53117d7d925SStefano Zampini    Level: intermediate
53217d7d925SStefano Zampini 
53317d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation()
53417d7d925SStefano Zampini @*/
535b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIInt(PetscInt n,PetscMPIInt X[])
53617d7d925SStefano Zampini {
53759e16d97SJunchao Zhang   PetscErrorCode ierr;
53859e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
53917d7d925SStefano Zampini 
54017d7d925SStefano Zampini   PetscFunctionBegin;
54159e16d97SJunchao Zhang   QuickSort1(PetscSortMPIInt,X,n,pivot,t1,ierr);
54217d7d925SStefano Zampini   PetscFunctionReturn(0);
54317d7d925SStefano Zampini }
54417d7d925SStefano Zampini 
54517d7d925SStefano Zampini /*@
54617d7d925SStefano Zampini    PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries
54717d7d925SStefano Zampini 
54817d7d925SStefano Zampini    Not Collective
54917d7d925SStefano Zampini 
55017d7d925SStefano Zampini    Input Parameters:
55117d7d925SStefano Zampini +  n  - number of values
55259e16d97SJunchao Zhang -  X  - array of integers
55317d7d925SStefano Zampini 
55417d7d925SStefano Zampini    Output Parameter:
55517d7d925SStefano Zampini .  n - number of non-redundant values
55617d7d925SStefano Zampini 
55717d7d925SStefano Zampini    Level: intermediate
55817d7d925SStefano Zampini 
55917d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
56017d7d925SStefano Zampini @*/
56159e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt X[])
56217d7d925SStefano Zampini {
56317d7d925SStefano Zampini   PetscErrorCode ierr;
56417d7d925SStefano Zampini   PetscInt       i,s = 0,N = *n, b = 0;
56517d7d925SStefano Zampini 
56617d7d925SStefano Zampini   PetscFunctionBegin;
56759e16d97SJunchao Zhang   ierr = PetscSortMPIInt(N,X);CHKERRQ(ierr);
56817d7d925SStefano Zampini   for (i=0; i<N-1; i++) {
56959e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
57059e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
571a297a907SKarl Rupp     } else s++;
57217d7d925SStefano Zampini   }
57317d7d925SStefano Zampini   *n = N - s;
57417d7d925SStefano Zampini   PetscFunctionReturn(0);
57517d7d925SStefano Zampini }
576c1f0200aSDmitry Karpeev 
5774d615ea0SBarry Smith /*@
5784d615ea0SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order;
5794d615ea0SBarry Smith        changes a second array to match the sorted first array.
5804d615ea0SBarry Smith 
5814d615ea0SBarry Smith    Not Collective
5824d615ea0SBarry Smith 
5834d615ea0SBarry Smith    Input Parameters:
5844d615ea0SBarry Smith +  n  - number of values
58559e16d97SJunchao Zhang .  X  - array of integers
58659e16d97SJunchao Zhang -  Y  - second array of integers
5874d615ea0SBarry Smith 
5884d615ea0SBarry Smith    Level: intermediate
5894d615ea0SBarry Smith 
5904d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
5914d615ea0SBarry Smith @*/
592b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt X[],PetscMPIInt Y[])
5934d615ea0SBarry Smith {
5944d615ea0SBarry Smith   PetscErrorCode ierr;
59559e16d97SJunchao Zhang   PetscMPIInt    pivot,t1,t2;
5964d615ea0SBarry Smith 
5974d615ea0SBarry Smith   PetscFunctionBegin;
59859e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithArray,X,Y,n,pivot,t1,t2,ierr);
5995569a785SJunchao Zhang   PetscFunctionReturn(0);
6005569a785SJunchao Zhang }
6015569a785SJunchao Zhang 
6025569a785SJunchao Zhang /*@
6035569a785SJunchao Zhang    PetscSortMPIIntWithIntArray - Sorts an array of MPI integers in place in increasing order;
6045569a785SJunchao Zhang        changes a second array of Petsc intergers to match the sorted first array.
6055569a785SJunchao Zhang 
6065569a785SJunchao Zhang    Not Collective
6075569a785SJunchao Zhang 
6085569a785SJunchao Zhang    Input Parameters:
6095569a785SJunchao Zhang +  n  - number of values
61059e16d97SJunchao Zhang .  X  - array of MPI integers
61159e16d97SJunchao Zhang -  Y  - second array of Petsc integers
6125569a785SJunchao Zhang 
6135569a785SJunchao Zhang    Level: intermediate
6145569a785SJunchao Zhang 
6155569a785SJunchao Zhang    Notes: this routine is useful when one needs to sort MPI ranks with other integer arrays.
6165569a785SJunchao Zhang 
6175569a785SJunchao Zhang .seealso: PetscSortMPIIntWithArray()
6185569a785SJunchao Zhang @*/
619b4f531b9SJunchao Zhang PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n,PetscMPIInt X[],PetscInt Y[])
6205569a785SJunchao Zhang {
6215569a785SJunchao Zhang   PetscErrorCode ierr;
62259e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
6235569a785SJunchao Zhang   PetscInt       t2;
6245569a785SJunchao Zhang 
6255569a785SJunchao Zhang   PetscFunctionBegin;
62659e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithIntArray,X,Y,n,pivot,t1,t2,ierr);
627e5c89e4eSSatish Balay   PetscFunctionReturn(0);
628e5c89e4eSSatish Balay }
629e5c89e4eSSatish Balay 
630e5c89e4eSSatish Balay /*@
631e5c89e4eSSatish Balay    PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
632e5c89e4eSSatish Balay        changes a second SCALAR array to match the sorted first INTEGER array.
633e5c89e4eSSatish Balay 
634e5c89e4eSSatish Balay    Not Collective
635e5c89e4eSSatish Balay 
636e5c89e4eSSatish Balay    Input Parameters:
637e5c89e4eSSatish Balay +  n  - number of values
63859e16d97SJunchao Zhang .  X  - array of integers
63959e16d97SJunchao Zhang -  Y  - second array of scalars
640e5c89e4eSSatish Balay 
641e5c89e4eSSatish Balay    Level: intermediate
642e5c89e4eSSatish Balay 
643e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
644e5c89e4eSSatish Balay @*/
645b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithScalarArray(PetscInt n,PetscInt X[],PetscScalar Y[])
646e5c89e4eSSatish Balay {
647e5c89e4eSSatish Balay   PetscErrorCode ierr;
64859e16d97SJunchao Zhang   PetscInt       pivot,t1;
64959e16d97SJunchao Zhang   PetscScalar    t2;
650e5c89e4eSSatish Balay 
651e5c89e4eSSatish Balay   PetscFunctionBegin;
65259e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithScalarArray,X,Y,n,pivot,t1,t2,ierr);
65317df18f2SToby Isaac   PetscFunctionReturn(0);
65417df18f2SToby Isaac }
65517df18f2SToby Isaac 
6566c2863d0SBarry Smith /*@C
65717df18f2SToby Isaac    PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order;
65817df18f2SToby Isaac        changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
65917df18f2SToby Isaac        provide workspace (the size of an element in the data array) to use when sorting.
66017df18f2SToby Isaac 
66117df18f2SToby Isaac    Not Collective
66217df18f2SToby Isaac 
66317df18f2SToby Isaac    Input Parameters:
66417df18f2SToby Isaac +  n  - number of values
66559e16d97SJunchao Zhang .  X  - array of integers
66659e16d97SJunchao Zhang .  Y  - second array of data
66717df18f2SToby Isaac .  size - sizeof elements in the data array in bytes
66859e16d97SJunchao Zhang -  t2   - workspace of "size" bytes used when sorting
66917df18f2SToby Isaac 
67017df18f2SToby Isaac    Level: intermediate
67117df18f2SToby Isaac 
67217df18f2SToby Isaac .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
67317df18f2SToby Isaac @*/
674b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithDataArray(PetscInt n,PetscInt X[],void *Y,size_t size,void *t2)
67517df18f2SToby Isaac {
67617df18f2SToby Isaac   PetscErrorCode ierr;
67759e16d97SJunchao Zhang   char           *YY = (char*)Y;
67859e16d97SJunchao Zhang   PetscInt       i,j,p,t1,pivot,hi=n-1,l,r;
67917df18f2SToby Isaac 
68017df18f2SToby Isaac   PetscFunctionBegin;
68117df18f2SToby Isaac   if (n<8) {
68259e16d97SJunchao Zhang     for (i=0; i<n; i++) {
68359e16d97SJunchao Zhang       pivot = X[i];
68459e16d97SJunchao Zhang       for (j=i+1; j<n; j++) {
68559e16d97SJunchao Zhang         if (pivot > X[j]) {
68659e16d97SJunchao Zhang           SWAP2Data(X[i],X[j],YY+size*i,YY+size*j,t1,t2,size);
68759e16d97SJunchao Zhang           pivot = X[i];
68817df18f2SToby Isaac         }
68917df18f2SToby Isaac       }
69017df18f2SToby Isaac     }
69117df18f2SToby Isaac   } else {
69259e16d97SJunchao Zhang     /* Two way partition */
69359e16d97SJunchao Zhang     p     = MEDIAN(X,hi);
69459e16d97SJunchao Zhang     pivot = X[p];
69559e16d97SJunchao Zhang     l     = 0;
69659e16d97SJunchao Zhang     r     = hi;
69759e16d97SJunchao Zhang     while(1) {
69859e16d97SJunchao Zhang       while (X[l] < pivot) l++;
69959e16d97SJunchao Zhang       while (X[r] > pivot) r--;
70059e16d97SJunchao Zhang       if (l >= r) {r++; break;}
70159e16d97SJunchao Zhang       SWAP2Data(X[l],X[r],YY+size*l,YY+size*r,t1,t2,size);
70259e16d97SJunchao Zhang       l++;
70359e16d97SJunchao Zhang       r--;
70459e16d97SJunchao Zhang     }
70559e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(l,X,Y,size,t2);CHKERRQ(ierr);
70659e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(hi-r+1,X+r,YY+size*r,size,t2);CHKERRQ(ierr);
70717df18f2SToby Isaac   }
70817df18f2SToby Isaac   PetscFunctionReturn(0);
70917df18f2SToby Isaac }
71017df18f2SToby Isaac 
71121e72a00SBarry Smith /*@
71221e72a00SBarry Smith    PetscMergeIntArray -     Merges two SORTED integer arrays, removes duplicate elements.
71321e72a00SBarry Smith 
71421e72a00SBarry Smith    Not Collective
71521e72a00SBarry Smith 
71621e72a00SBarry Smith    Input Parameters:
71721e72a00SBarry Smith +  an  - number of values in the first array
71821e72a00SBarry Smith .  aI  - first sorted array of integers
71921e72a00SBarry Smith .  bn  - number of values in the second array
72021e72a00SBarry Smith -  bI  - second array of integers
72121e72a00SBarry Smith 
72221e72a00SBarry Smith    Output Parameters:
72321e72a00SBarry Smith +  n   - number of values in the merged array
7246c2863d0SBarry Smith -  L   - merged sorted array, this is allocated if an array is not provided
72521e72a00SBarry Smith 
72621e72a00SBarry Smith    Level: intermediate
72721e72a00SBarry Smith 
72821e72a00SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
72921e72a00SBarry Smith @*/
7306c2863d0SBarry Smith PetscErrorCode  PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
73121e72a00SBarry Smith {
73221e72a00SBarry Smith   PetscErrorCode ierr;
73321e72a00SBarry Smith   PetscInt       *L_ = *L, ak, bk, k;
73421e72a00SBarry Smith 
73521e72a00SBarry Smith   if (!L_) {
73621e72a00SBarry Smith     ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr);
73721e72a00SBarry Smith     L_   = *L;
73821e72a00SBarry Smith   }
73921e72a00SBarry Smith   k = ak = bk = 0;
74021e72a00SBarry Smith   while (ak < an && bk < bn) {
74121e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
74221e72a00SBarry Smith       L_[k] = aI[ak];
74321e72a00SBarry Smith       ++ak;
74421e72a00SBarry Smith       ++bk;
74521e72a00SBarry Smith       ++k;
74621e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
74721e72a00SBarry Smith       L_[k] = aI[ak];
74821e72a00SBarry Smith       ++ak;
74921e72a00SBarry Smith       ++k;
75021e72a00SBarry Smith     } else {
75121e72a00SBarry Smith       L_[k] = bI[bk];
75221e72a00SBarry Smith       ++bk;
75321e72a00SBarry Smith       ++k;
75421e72a00SBarry Smith     }
75521e72a00SBarry Smith   }
75621e72a00SBarry Smith   if (ak < an) {
757580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
75821e72a00SBarry Smith     k   += (an-ak);
75921e72a00SBarry Smith   }
76021e72a00SBarry Smith   if (bk < bn) {
761580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
76221e72a00SBarry Smith     k   += (bn-bk);
76321e72a00SBarry Smith   }
76421e72a00SBarry Smith   *n = k;
76521e72a00SBarry Smith   PetscFunctionReturn(0);
76621e72a00SBarry Smith }
767b4301105SBarry Smith 
768c1f0200aSDmitry Karpeev /*@
76921e72a00SBarry Smith    PetscMergeIntArrayPair -     Merges two SORTED integer arrays that share NO common values along with an additional array of integers.
770c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
771c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
772c1f0200aSDmitry Karpeev 
773c1f0200aSDmitry Karpeev    Not Collective
774c1f0200aSDmitry Karpeev 
775c1f0200aSDmitry Karpeev    Input Parameters:
776c1f0200aSDmitry Karpeev +  an  - number of values in the first array
777c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
778c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
779c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
780c1f0200aSDmitry Karpeev .  bI  - second array of integers
781c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
782c1f0200aSDmitry Karpeev 
783c1f0200aSDmitry Karpeev    Output Parameters:
784c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
78514c49006SJed Brown .  L   - merged sorted array
786c1f0200aSDmitry Karpeev -  J   - merged additional array
787c1f0200aSDmitry Karpeev 
78895452b02SPatrick Sanan    Notes:
78995452b02SPatrick 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
790c1f0200aSDmitry Karpeev    Level: intermediate
791c1f0200aSDmitry Karpeev 
792c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
793c1f0200aSDmitry Karpeev @*/
7946c2863d0SBarry 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)
795c1f0200aSDmitry Karpeev {
796c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
79770d8d27cSBarry Smith   PetscInt       n_, *L_, *J_, ak, bk, k;
798c1f0200aSDmitry Karpeev 
79914c49006SJed Brown   PetscFunctionBegin;
80070d8d27cSBarry Smith   PetscValidIntPointer(L,8);
80170d8d27cSBarry Smith   PetscValidIntPointer(J,9);
802c1f0200aSDmitry Karpeev   n_ = an + bn;
803c1f0200aSDmitry Karpeev   *n = n_;
80470d8d27cSBarry Smith   if (!*L) {
805785e854fSJed Brown     ierr = PetscMalloc1(n_, L);CHKERRQ(ierr);
80670d8d27cSBarry Smith   }
807d7aa01a8SHong Zhang   L_ = *L;
80870d8d27cSBarry Smith   if (!*J) {
80970d8d27cSBarry Smith     ierr = PetscMalloc1(n_, J);CHKERRQ(ierr);
810c1f0200aSDmitry Karpeev   }
811c1f0200aSDmitry Karpeev   J_   = *J;
812c1f0200aSDmitry Karpeev   k = ak = bk = 0;
813c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
814c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
815d7aa01a8SHong Zhang       L_[k] = aI[ak];
816c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
817c1f0200aSDmitry Karpeev       ++ak;
818c1f0200aSDmitry Karpeev       ++k;
819a297a907SKarl Rupp     } else {
820d7aa01a8SHong Zhang       L_[k] = bI[bk];
821c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
822c1f0200aSDmitry Karpeev       ++bk;
823c1f0200aSDmitry Karpeev       ++k;
824c1f0200aSDmitry Karpeev     }
825c1f0200aSDmitry Karpeev   }
826c1f0200aSDmitry Karpeev   if (ak < an) {
827580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
828580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,aJ+ak,an-ak);CHKERRQ(ierr);
829c1f0200aSDmitry Karpeev     k   += (an-ak);
830c1f0200aSDmitry Karpeev   }
831c1f0200aSDmitry Karpeev   if (bk < bn) {
832580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
833580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,bJ+bk,bn-bk);CHKERRQ(ierr);
834c1f0200aSDmitry Karpeev   }
835c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
836c1f0200aSDmitry Karpeev }
837c1f0200aSDmitry Karpeev 
838e498c390SJed Brown /*@
839e498c390SJed Brown    PetscMergeMPIIntArray -     Merges two SORTED integer arrays.
840e498c390SJed Brown 
841e498c390SJed Brown    Not Collective
842e498c390SJed Brown 
843e498c390SJed Brown    Input Parameters:
844e498c390SJed Brown +  an  - number of values in the first array
845e498c390SJed Brown .  aI  - first sorted array of integers
846e498c390SJed Brown .  bn  - number of values in the second array
847e498c390SJed Brown -  bI  - second array of integers
848e498c390SJed Brown 
849e498c390SJed Brown    Output Parameters:
850e498c390SJed Brown +  n   - number of values in the merged array (<= an + bn)
851e498c390SJed Brown -  L   - merged sorted array, allocated if address of NULL pointer is passed
852e498c390SJed Brown 
853e498c390SJed Brown    Level: intermediate
854e498c390SJed Brown 
855e498c390SJed Brown .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
856e498c390SJed Brown @*/
857e498c390SJed Brown PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L)
858e498c390SJed Brown {
859e498c390SJed Brown   PetscErrorCode ierr;
860e498c390SJed Brown   PetscInt       ai,bi,k;
861e498c390SJed Brown 
862e498c390SJed Brown   PetscFunctionBegin;
863e498c390SJed Brown   if (!*L) {ierr = PetscMalloc1((an+bn),L);CHKERRQ(ierr);}
864e498c390SJed Brown   for (ai=0,bi=0,k=0; ai<an || bi<bn; ) {
865e498c390SJed Brown     PetscInt t = -1;
866e498c390SJed Brown     for ( ; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
867e498c390SJed Brown     for ( ; bi<bn && bI[bi] == t; bi++);
868e498c390SJed Brown     for ( ; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
869e498c390SJed Brown     for ( ; ai<an && aI[ai] == t; ai++);
870e498c390SJed Brown   }
871e498c390SJed Brown   *n = k;
872e498c390SJed Brown   PetscFunctionReturn(0);
873e498c390SJed Brown }
874e5c89e4eSSatish Balay 
8756c2863d0SBarry Smith /*@C
87642eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
87742eaa7f3SBarry Smith 
87842eaa7f3SBarry Smith    Not Collective
87942eaa7f3SBarry Smith 
88042eaa7f3SBarry Smith    Input Parameters:
88142eaa7f3SBarry Smith +  n  - number of values
88242eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
88342eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
88442eaa7f3SBarry Smith 
88542eaa7f3SBarry Smith    Output Parameters:
88642eaa7f3SBarry Smith +  Nlevels - the number of levels
88742eaa7f3SBarry Smith .  Level - for each node tells its level
88842eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
88942eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
89042eaa7f3SBarry Smith -  Column - for each id tells its column index
89142eaa7f3SBarry Smith 
8926c2863d0SBarry Smith    Level: developer
89342eaa7f3SBarry Smith 
89495452b02SPatrick Sanan    Notes:
89595452b02SPatrick Sanan     This code is not currently used
89642eaa7f3SBarry Smith 
89742eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation()
89842eaa7f3SBarry Smith @*/
8997087cfbeSBarry Smith PetscErrorCode  PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column)
90042eaa7f3SBarry Smith {
90142eaa7f3SBarry Smith   PetscInt       i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column;
90242eaa7f3SBarry Smith   PetscErrorCode ierr;
903ace3abfcSBarry Smith   PetscBool      done = PETSC_FALSE;
90442eaa7f3SBarry Smith 
90542eaa7f3SBarry Smith   PetscFunctionBegin;
90642eaa7f3SBarry Smith   if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set");
90742eaa7f3SBarry Smith   for (i=0; i<n; i++) {
90842eaa7f3SBarry Smith     if (mask[i]) continue;
90942eaa7f3SBarry Smith     if (parentid[i]  == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent");
91042eaa7f3SBarry Smith     if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked");
91142eaa7f3SBarry Smith   }
91242eaa7f3SBarry Smith 
91342eaa7f3SBarry Smith   for (i=0; i<n; i++) {
91442eaa7f3SBarry Smith     if (!mask[i]) nmask++;
91542eaa7f3SBarry Smith   }
91642eaa7f3SBarry Smith 
91742eaa7f3SBarry Smith   /* determine the level in the tree of each node */
9181795a4d1SJed Brown   ierr = PetscCalloc1(n,&level);CHKERRQ(ierr);
919a297a907SKarl Rupp 
92042eaa7f3SBarry Smith   level[0] = 1;
92142eaa7f3SBarry Smith   while (!done) {
92242eaa7f3SBarry Smith     done = PETSC_TRUE;
92342eaa7f3SBarry Smith     for (i=0; i<n; i++) {
92442eaa7f3SBarry Smith       if (mask[i]) continue;
92542eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
92642eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
92742eaa7f3SBarry Smith     }
92842eaa7f3SBarry Smith   }
92942eaa7f3SBarry Smith   for (i=0; i<n; i++) {
93042eaa7f3SBarry Smith     level[i]--;
93142eaa7f3SBarry Smith     nlevels = PetscMax(nlevels,level[i]);
93242eaa7f3SBarry Smith   }
93342eaa7f3SBarry Smith 
93442eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
9351795a4d1SJed Brown   ierr = PetscCalloc1(nlevels,&levelcnt);CHKERRQ(ierr);
93642eaa7f3SBarry Smith   for (i=0; i<n; i++) {
93742eaa7f3SBarry Smith     if (mask[i]) continue;
93842eaa7f3SBarry Smith     levelcnt[level[i]-1]++;
93942eaa7f3SBarry Smith   }
940a297a907SKarl Rupp   for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]);
94142eaa7f3SBarry Smith 
94242eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
943dcca6d9dSJed Brown   ierr = PetscMalloc2(levelmax,&workid,levelmax,&workparentid);CHKERRQ(ierr);
944785e854fSJed Brown   ierr = PetscMalloc1(nmask,&idbylevel);CHKERRQ(ierr);
94542eaa7f3SBarry Smith   for (j=1; j<=nlevels;j++) {
94642eaa7f3SBarry Smith     cnt = 0;
94742eaa7f3SBarry Smith     for (i=0; i<n; i++) {
94842eaa7f3SBarry Smith       if (mask[i]) continue;
94942eaa7f3SBarry Smith       if (level[i] != j) continue;
95042eaa7f3SBarry Smith       workid[cnt]         = i;
95142eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
95242eaa7f3SBarry Smith     }
95342eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
95442eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
95542eaa7f3SBarry Smith     ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr);
95642eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
95742eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
958580bdb30SBarry Smith     ierr  = PetscArraycpy(idbylevel+tcnt,workid,cnt);CHKERRQ(ierr);
95942eaa7f3SBarry Smith     tcnt += cnt;
96042eaa7f3SBarry Smith   }
96142eaa7f3SBarry Smith   if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes");
96242eaa7f3SBarry Smith   ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr);
96342eaa7f3SBarry Smith 
96442eaa7f3SBarry Smith   /* for each node list its column */
965785e854fSJed Brown   ierr = PetscMalloc1(n,&column);CHKERRQ(ierr);
96642eaa7f3SBarry Smith   cnt = 0;
96742eaa7f3SBarry Smith   for (j=0; j<nlevels; j++) {
96842eaa7f3SBarry Smith     for (i=0; i<levelcnt[j]; i++) {
96942eaa7f3SBarry Smith       column[idbylevel[cnt++]] = i;
97042eaa7f3SBarry Smith     }
97142eaa7f3SBarry Smith   }
97242eaa7f3SBarry Smith 
97342eaa7f3SBarry Smith   *Nlevels   = nlevels;
97442eaa7f3SBarry Smith   *Level     = level;
97542eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
97642eaa7f3SBarry Smith   *Idbylevel = idbylevel;
97742eaa7f3SBarry Smith   *Column    = column;
97842eaa7f3SBarry Smith   PetscFunctionReturn(0);
97942eaa7f3SBarry Smith }
980ce605777SToby Isaac 
981ce605777SToby Isaac /*@
982ce605777SToby Isaac   PetscParallelSortedInt - Check whether an integer array, distributed over a communicator, is globally sorted.
983ce605777SToby Isaac 
984ce605777SToby Isaac   Collective
985ce605777SToby Isaac 
986ce605777SToby Isaac   Input Parameters:
987ce605777SToby Isaac + comm - the MPI communicator
988ce605777SToby Isaac . n - the local number of integers
989ce605777SToby Isaac - keys - the local array of integers
990ce605777SToby Isaac 
991ce605777SToby Isaac   Output Parameters:
992ce605777SToby Isaac . is_sorted - whether the array is globally sorted
993ce605777SToby Isaac 
994ce605777SToby Isaac   Level: developer
995ce605777SToby Isaac 
996ce605777SToby Isaac .seealso: PetscParallelSortInt()
997ce605777SToby Isaac @*/
998ce605777SToby Isaac PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
999ce605777SToby Isaac {
1000ce605777SToby Isaac   PetscBool      sorted;
1001ce605777SToby Isaac   PetscInt       i, min, max, prevmax;
10022a1da528SToby Isaac   PetscMPIInt    rank;
1003ce605777SToby Isaac   PetscErrorCode ierr;
1004ce605777SToby Isaac 
1005ce605777SToby Isaac   PetscFunctionBegin;
1006ce605777SToby Isaac   sorted = PETSC_TRUE;
1007ce605777SToby Isaac   min    = PETSC_MAX_INT;
1008ce605777SToby Isaac   max    = PETSC_MIN_INT;
1009ce605777SToby Isaac   if (n) {
1010ce605777SToby Isaac     min = keys[0];
1011ce605777SToby Isaac     max = keys[0];
1012ce605777SToby Isaac   }
1013ce605777SToby Isaac   for (i = 1; i < n; i++) {
1014ce605777SToby Isaac     if (keys[i] < keys[i - 1]) break;
1015ce605777SToby Isaac     min = PetscMin(min,keys[i]);
1016ce605777SToby Isaac     max = PetscMax(max,keys[i]);
1017ce605777SToby Isaac   }
1018ce605777SToby Isaac   if (i < n) sorted = PETSC_FALSE;
1019ce605777SToby Isaac   prevmax = PETSC_MIN_INT;
1020ce605777SToby Isaac   ierr = MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
10212a1da528SToby Isaac   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10222a1da528SToby Isaac   if (!rank) prevmax = PETSC_MIN_INT;
1023ce605777SToby Isaac   if (prevmax > min) sorted = PETSC_FALSE;
1024ce605777SToby Isaac   ierr = MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr);
1025ce605777SToby Isaac   PetscFunctionReturn(0);
1026ce605777SToby Isaac }
1027ce605777SToby Isaac 
1028