xref: /petsc/src/sys/utils/sorti.c (revision 981bb840e9420dc4902ca8e1ce5c98c517873cd2)
17d0a6c19SBarry Smith 
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file contains routines for sorting integers. Values are sorted in place.
4c4762a1bSJed Brown    One can use src/sys/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 {                                                                           \
126*981bb840SJunchao Zhang     PetscCount 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 {                                                                           \
149*981bb840SJunchao Zhang     PetscCount 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 {                                                                           \
171*981bb840SJunchao Zhang     PetscCount 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 {                                                                           \
193*981bb840SJunchao Zhang     PetscCount 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 
245676f2a66SJacob Faibussowitsch    Notes:
246676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscIntSortSemiOrdered(), and may perform faster especially if the array
247a5b23f4aSJose E. Roman    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
248676f2a66SJacob Faibussowitsch    code to see which routine is fastest.
249676f2a66SJacob Faibussowitsch 
250e5c89e4eSSatish Balay    Level: intermediate
251e5c89e4eSSatish Balay 
252676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation()
253e5c89e4eSSatish Balay @*/
254b4f531b9SJunchao Zhang PetscErrorCode  PetscSortInt(PetscInt n,PetscInt X[])
255e5c89e4eSSatish Balay {
25659e16d97SJunchao Zhang   PetscErrorCode ierr;
25759e16d97SJunchao Zhang   PetscInt       pivot,t1;
258e5c89e4eSSatish Balay 
259e5c89e4eSSatish Balay   PetscFunctionBegin;
26059e16d97SJunchao Zhang   QuickSort1(PetscSortInt,X,n,pivot,t1,ierr);
261e5c89e4eSSatish Balay   PetscFunctionReturn(0);
262e5c89e4eSSatish Balay }
263e5c89e4eSSatish Balay 
2641db36a52SBarry Smith /*@
265ce605777SToby Isaac    PetscSortReverseInt - Sorts an array of integers in place in decreasing order.
266ce605777SToby Isaac 
267ce605777SToby Isaac    Not Collective
268ce605777SToby Isaac 
269ce605777SToby Isaac    Input Parameters:
270ce605777SToby Isaac +  n  - number of values
271ce605777SToby Isaac -  X  - array of integers
272ce605777SToby Isaac 
273ce605777SToby Isaac    Level: intermediate
274ce605777SToby Isaac 
275676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithPermutation()
276ce605777SToby Isaac @*/
277ce605777SToby Isaac PetscErrorCode  PetscSortReverseInt(PetscInt n,PetscInt X[])
278ce605777SToby Isaac {
279ce605777SToby Isaac   PetscErrorCode ierr;
280ce605777SToby Isaac   PetscInt       pivot,t1;
281ce605777SToby Isaac 
282ce605777SToby Isaac   PetscFunctionBegin;
283ce605777SToby Isaac   QuickSortReverse1(PetscSortReverseInt,X,n,pivot,t1,ierr);
284ce605777SToby Isaac   PetscFunctionReturn(0);
285ce605777SToby Isaac }
286ce605777SToby Isaac 
287ce605777SToby Isaac /*@
28822ab5688SLisandro Dalcin    PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted input array
28922ab5688SLisandro Dalcin 
29022ab5688SLisandro Dalcin    Not Collective
29122ab5688SLisandro Dalcin 
29222ab5688SLisandro Dalcin    Input Parameters:
29322ab5688SLisandro Dalcin +  n  - number of values
29459e16d97SJunchao Zhang -  X  - sorted array of integers
29522ab5688SLisandro Dalcin 
29622ab5688SLisandro Dalcin    Output Parameter:
29722ab5688SLisandro Dalcin .  n - number of non-redundant values
29822ab5688SLisandro Dalcin 
29922ab5688SLisandro Dalcin    Level: intermediate
30022ab5688SLisandro Dalcin 
30122ab5688SLisandro Dalcin .seealso: PetscSortInt()
30222ab5688SLisandro Dalcin @*/
30359e16d97SJunchao Zhang PetscErrorCode  PetscSortedRemoveDupsInt(PetscInt *n,PetscInt X[])
30422ab5688SLisandro Dalcin {
30522ab5688SLisandro Dalcin   PetscInt i,s = 0,N = *n, b = 0;
30622ab5688SLisandro Dalcin 
30722ab5688SLisandro Dalcin   PetscFunctionBegin;
308fb61b9e4SStefano Zampini   PetscCheckSorted(*n,X);
30922ab5688SLisandro Dalcin   for (i=0; i<N-1; i++) {
31059e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
31159e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
31222ab5688SLisandro Dalcin     } else s++;
31322ab5688SLisandro Dalcin   }
31422ab5688SLisandro Dalcin   *n = N - s;
31522ab5688SLisandro Dalcin   PetscFunctionReturn(0);
31622ab5688SLisandro Dalcin }
31722ab5688SLisandro Dalcin 
31822ab5688SLisandro Dalcin /*@
3191db36a52SBarry Smith    PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries
3201db36a52SBarry Smith 
3211db36a52SBarry Smith    Not Collective
3221db36a52SBarry Smith 
3231db36a52SBarry Smith    Input Parameters:
3241db36a52SBarry Smith +  n  - number of values
32559e16d97SJunchao Zhang -  X  - array of integers
3261db36a52SBarry Smith 
3271db36a52SBarry Smith    Output Parameter:
3281db36a52SBarry Smith .  n - number of non-redundant values
3291db36a52SBarry Smith 
3301db36a52SBarry Smith    Level: intermediate
3311db36a52SBarry Smith 
332676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortedRemoveDupsInt()
3331db36a52SBarry Smith @*/
33459e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsInt(PetscInt *n,PetscInt X[])
3351db36a52SBarry Smith {
3361db36a52SBarry Smith   PetscErrorCode ierr;
3371db36a52SBarry Smith 
3381db36a52SBarry Smith   PetscFunctionBegin;
33959e16d97SJunchao Zhang   ierr = PetscSortInt(*n,X);CHKERRQ(ierr);
34059e16d97SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(n,X);CHKERRQ(ierr);
3411db36a52SBarry Smith   PetscFunctionReturn(0);
3421db36a52SBarry Smith }
3431db36a52SBarry Smith 
34460e03357SMatthew G Knepley /*@
345d9f1162dSJed Brown   PetscFindInt - Finds integer in a sorted array of integers
34660e03357SMatthew G Knepley 
34760e03357SMatthew G Knepley    Not Collective
34860e03357SMatthew G Knepley 
34960e03357SMatthew G Knepley    Input Parameters:
35060e03357SMatthew G Knepley +  key - the integer to locate
351d9f1162dSJed Brown .  n   - number of values in the array
35259e16d97SJunchao Zhang -  X  - array of integers
35360e03357SMatthew G Knepley 
35460e03357SMatthew G Knepley    Output Parameter:
355d9f1162dSJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
35660e03357SMatthew G Knepley 
35760e03357SMatthew G Knepley    Level: intermediate
35860e03357SMatthew G Knepley 
359676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
36060e03357SMatthew G Knepley @*/
36159e16d97SJunchao Zhang PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
36260e03357SMatthew G Knepley {
363d9f1162dSJed Brown   PetscInt lo = 0,hi = n;
36460e03357SMatthew G Knepley 
36560e03357SMatthew G Knepley   PetscFunctionBegin;
36660e03357SMatthew G Knepley   PetscValidPointer(loc,4);
36798781d82SMatthew G Knepley   if (!n) {*loc = -1; PetscFunctionReturn(0);}
36859e16d97SJunchao Zhang   PetscValidPointer(X,3);
3696a7c706bSVaclav Hapla   PetscCheckSorted(n,X);
370d9f1162dSJed Brown   while (hi - lo > 1) {
371d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo)/2;
37259e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
373d9f1162dSJed Brown     else               lo = mid;
37460e03357SMatthew G Knepley   }
37559e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
37660e03357SMatthew G Knepley   PetscFunctionReturn(0);
37760e03357SMatthew G Knepley }
37860e03357SMatthew G Knepley 
379d2aeb606SJed Brown /*@
380f1cab4e1SJunchao Zhang   PetscCheckDupsInt - Checks if an integer array has duplicates
381f1cab4e1SJunchao Zhang 
382f1cab4e1SJunchao Zhang    Not Collective
383f1cab4e1SJunchao Zhang 
384f1cab4e1SJunchao Zhang    Input Parameters:
385f1cab4e1SJunchao Zhang +  n  - number of values in the array
386f1cab4e1SJunchao Zhang -  X  - array of integers
387f1cab4e1SJunchao Zhang 
388f1cab4e1SJunchao Zhang    Output Parameter:
389f1cab4e1SJunchao Zhang .  dups - True if the array has dups, otherwise false
390f1cab4e1SJunchao Zhang 
391f1cab4e1SJunchao Zhang    Level: intermediate
392f1cab4e1SJunchao Zhang 
393f1cab4e1SJunchao Zhang .seealso: PetscSortRemoveDupsInt()
394f1cab4e1SJunchao Zhang @*/
395f1cab4e1SJunchao Zhang PetscErrorCode PetscCheckDupsInt(PetscInt n,const PetscInt X[],PetscBool *dups)
396f1cab4e1SJunchao Zhang {
397f1cab4e1SJunchao Zhang   PetscErrorCode ierr;
398f1cab4e1SJunchao Zhang   PetscInt       i;
399f1cab4e1SJunchao Zhang   PetscHSetI     ht;
400f1cab4e1SJunchao Zhang   PetscBool      missing;
401f1cab4e1SJunchao Zhang 
402f1cab4e1SJunchao Zhang   PetscFunctionBegin;
403f1cab4e1SJunchao Zhang   PetscValidPointer(dups,3);
404f1cab4e1SJunchao Zhang   *dups = PETSC_FALSE;
405f1cab4e1SJunchao Zhang   if (n > 1) {
406f1cab4e1SJunchao Zhang     ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
407f1cab4e1SJunchao Zhang     ierr = PetscHSetIResize(ht,n);CHKERRQ(ierr);
408f1cab4e1SJunchao Zhang     for (i=0; i<n; i++) {
409f1cab4e1SJunchao Zhang       ierr = PetscHSetIQueryAdd(ht,X[i],&missing);CHKERRQ(ierr);
410f1cab4e1SJunchao Zhang       if (!missing) {*dups = PETSC_TRUE; break;}
411f1cab4e1SJunchao Zhang     }
412f1cab4e1SJunchao Zhang     ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
413f1cab4e1SJunchao Zhang   }
414f1cab4e1SJunchao Zhang   PetscFunctionReturn(0);
415f1cab4e1SJunchao Zhang }
416f1cab4e1SJunchao Zhang 
417f1cab4e1SJunchao Zhang /*@
418d2aeb606SJed Brown   PetscFindMPIInt - Finds MPI integer in a sorted array of integers
419d2aeb606SJed Brown 
420d2aeb606SJed Brown    Not Collective
421d2aeb606SJed Brown 
422d2aeb606SJed Brown    Input Parameters:
423d2aeb606SJed Brown +  key - the integer to locate
424d2aeb606SJed Brown .  n   - number of values in the array
42559e16d97SJunchao Zhang -  X   - array of integers
426d2aeb606SJed Brown 
427d2aeb606SJed Brown    Output Parameter:
428d2aeb606SJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
429d2aeb606SJed Brown 
430d2aeb606SJed Brown    Level: intermediate
431d2aeb606SJed Brown 
432676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
433d2aeb606SJed Brown @*/
43459e16d97SJunchao Zhang PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
435d2aeb606SJed Brown {
436d2aeb606SJed Brown   PetscInt lo = 0,hi = n;
437d2aeb606SJed Brown 
438d2aeb606SJed Brown   PetscFunctionBegin;
439d2aeb606SJed Brown   PetscValidPointer(loc,4);
440d2aeb606SJed Brown   if (!n) {*loc = -1; PetscFunctionReturn(0);}
44159e16d97SJunchao Zhang   PetscValidPointer(X,3);
4426a7c706bSVaclav Hapla   PetscCheckSorted(n,X);
443d2aeb606SJed Brown   while (hi - lo > 1) {
444d2aeb606SJed Brown     PetscInt mid = lo + (hi - lo)/2;
44559e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
446d2aeb606SJed Brown     else               lo = mid;
447d2aeb606SJed Brown   }
44859e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
449e5c89e4eSSatish Balay   PetscFunctionReturn(0);
450e5c89e4eSSatish Balay }
451e5c89e4eSSatish Balay 
452e5c89e4eSSatish Balay /*@
453e5c89e4eSSatish Balay    PetscSortIntWithArray - Sorts an array of integers in place in increasing order;
454e5c89e4eSSatish Balay        changes a second array to match the sorted first array.
455e5c89e4eSSatish Balay 
456e5c89e4eSSatish Balay    Not Collective
457e5c89e4eSSatish Balay 
458e5c89e4eSSatish Balay    Input Parameters:
459e5c89e4eSSatish Balay +  n  - number of values
46059e16d97SJunchao Zhang .  X  - array of integers
46159e16d97SJunchao Zhang -  Y  - second array of integers
462e5c89e4eSSatish Balay 
463e5c89e4eSSatish Balay    Level: intermediate
464e5c89e4eSSatish Balay 
465*981bb840SJunchao Zhang .seealso: PetscIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithCountArray()
466e5c89e4eSSatish Balay @*/
467b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithArray(PetscInt n,PetscInt X[],PetscInt Y[])
468e5c89e4eSSatish Balay {
469e5c89e4eSSatish Balay   PetscErrorCode ierr;
47059e16d97SJunchao Zhang   PetscInt       pivot,t1,t2;
471e5c89e4eSSatish Balay 
472e5c89e4eSSatish Balay   PetscFunctionBegin;
47359e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithArray,X,Y,n,pivot,t1,t2,ierr);
474c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
475c1f0200aSDmitry Karpeev }
476c1f0200aSDmitry Karpeev 
477c1f0200aSDmitry Karpeev /*@
478c1f0200aSDmitry Karpeev    PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order;
479c1f0200aSDmitry Karpeev        changes a pair of integer arrays to match the sorted first array.
480c1f0200aSDmitry Karpeev 
481c1f0200aSDmitry Karpeev    Not Collective
482c1f0200aSDmitry Karpeev 
483c1f0200aSDmitry Karpeev    Input Parameters:
484c1f0200aSDmitry Karpeev +  n  - number of values
48559e16d97SJunchao Zhang .  X  - array of integers
48659e16d97SJunchao Zhang .  Y  - second array of integers (first array of the pair)
48759e16d97SJunchao Zhang -  Z  - third array of integers  (second array of the pair)
488c1f0200aSDmitry Karpeev 
489c1f0200aSDmitry Karpeev    Level: intermediate
490c1f0200aSDmitry Karpeev 
491*981bb840SJunchao Zhang .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortIntWithArray(), PetscIntSortSemiOrdered(), PetscSortIntWithIntCountArrayPair()
492c1f0200aSDmitry Karpeev @*/
493b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithArrayPair(PetscInt n,PetscInt X[],PetscInt Y[],PetscInt Z[])
494c1f0200aSDmitry Karpeev {
495c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
49659e16d97SJunchao Zhang   PetscInt       pivot,t1,t2,t3;
497c1f0200aSDmitry Karpeev 
498c1f0200aSDmitry Karpeev   PetscFunctionBegin;
49959e16d97SJunchao Zhang   QuickSort3(PetscSortIntWithArrayPair,X,Y,Z,n,pivot,t1,t2,t3,ierr);
500c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
501c1f0200aSDmitry Karpeev }
502c1f0200aSDmitry Karpeev 
50317d7d925SStefano Zampini /*@
504*981bb840SJunchao Zhang    PetscSortIntWithCountArray - Sorts an array of integers in place in increasing order;
505*981bb840SJunchao Zhang        changes a second array to match the sorted first array.
506*981bb840SJunchao Zhang 
507*981bb840SJunchao Zhang    Not Collective
508*981bb840SJunchao Zhang 
509*981bb840SJunchao Zhang    Input Parameters:
510*981bb840SJunchao Zhang +  n  - number of values
511*981bb840SJunchao Zhang .  X  - array of integers
512*981bb840SJunchao Zhang -  Y  - second array of PetscCounts (signed integers)
513*981bb840SJunchao Zhang 
514*981bb840SJunchao Zhang    Level: intermediate
515*981bb840SJunchao Zhang 
516*981bb840SJunchao Zhang .seealso: PetscIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
517*981bb840SJunchao Zhang @*/
518*981bb840SJunchao Zhang PetscErrorCode  PetscSortIntWithCountArray(PetscCount n,PetscInt X[],PetscCount Y[])
519*981bb840SJunchao Zhang {
520*981bb840SJunchao Zhang   PetscErrorCode ierr;
521*981bb840SJunchao Zhang   PetscInt       pivot,t1;
522*981bb840SJunchao Zhang   PetscCount     t2;
523*981bb840SJunchao Zhang 
524*981bb840SJunchao Zhang   PetscFunctionBegin;
525*981bb840SJunchao Zhang   QuickSort2(PetscSortIntWithCountArray,X,Y,n,pivot,t1,t2,ierr);
526*981bb840SJunchao Zhang   PetscFunctionReturn(0);
527*981bb840SJunchao Zhang }
528*981bb840SJunchao Zhang 
529*981bb840SJunchao Zhang /*@
530*981bb840SJunchao Zhang    PetscSortIntWithIntCountArrayPair - Sorts an array of integers in place in increasing order;
531*981bb840SJunchao Zhang        changes an integer array and a PetscCount array to match the sorted first array.
532*981bb840SJunchao Zhang 
533*981bb840SJunchao Zhang    Not Collective
534*981bb840SJunchao Zhang 
535*981bb840SJunchao Zhang    Input Parameters:
536*981bb840SJunchao Zhang +  n  - number of values
537*981bb840SJunchao Zhang .  X  - array of integers
538*981bb840SJunchao Zhang .  Y  - second array of integers (first array of the pair)
539*981bb840SJunchao Zhang -  Z  - third array of PetscCounts  (second array of the pair)
540*981bb840SJunchao Zhang 
541*981bb840SJunchao Zhang    Level: intermediate
542*981bb840SJunchao Zhang 
543*981bb840SJunchao Zhang    Notes:
544*981bb840SJunchao Zhang     Usually X, Y are matrix row/column indices, and Z is a permutation array and therefore Z's type is PetscCount to allow 2B+ nonzeros even with 32-bit PetscInt.
545*981bb840SJunchao Zhang 
546*981bb840SJunchao Zhang .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray(), PetscIntSortSemiOrdered(), PetscSortIntWithArrayPair()
547*981bb840SJunchao Zhang @*/
548*981bb840SJunchao Zhang PetscErrorCode  PetscSortIntWithIntCountArrayPair(PetscCount n,PetscInt X[],PetscInt Y[],PetscCount Z[])
549*981bb840SJunchao Zhang {
550*981bb840SJunchao Zhang   PetscErrorCode ierr;
551*981bb840SJunchao Zhang   PetscInt       pivot,t1,t2; /* pivot is take from X[], so its type is still PetscInt */
552*981bb840SJunchao Zhang   PetscCount     t3; /* temp for Z[] */
553*981bb840SJunchao Zhang 
554*981bb840SJunchao Zhang   PetscFunctionBegin;
555*981bb840SJunchao Zhang   QuickSort3(PetscSortIntWithIntCountArrayPair,X,Y,Z,n,pivot,t1,t2,t3,ierr);
556*981bb840SJunchao Zhang   PetscFunctionReturn(0);
557*981bb840SJunchao Zhang }
558*981bb840SJunchao Zhang 
559*981bb840SJunchao Zhang /*@
5606a7c706bSVaclav Hapla    PetscSortedMPIInt - Determines whether the array is sorted.
5616a7c706bSVaclav Hapla 
5626a7c706bSVaclav Hapla    Not Collective
5636a7c706bSVaclav Hapla 
5646a7c706bSVaclav Hapla    Input Parameters:
5656a7c706bSVaclav Hapla +  n  - number of values
5666a7c706bSVaclav Hapla -  X  - array of integers
5676a7c706bSVaclav Hapla 
5686a7c706bSVaclav Hapla    Output Parameters:
5696a7c706bSVaclav Hapla .  sorted - flag whether the array is sorted
5706a7c706bSVaclav Hapla 
5716a7c706bSVaclav Hapla    Level: intermediate
5726a7c706bSVaclav Hapla 
573676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortMPIInt(), PetscSortedInt(), PetscSortedReal()
5746a7c706bSVaclav Hapla @*/
5756a7c706bSVaclav Hapla PetscErrorCode  PetscSortedMPIInt(PetscInt n,const PetscMPIInt X[],PetscBool *sorted)
5766a7c706bSVaclav Hapla {
5776a7c706bSVaclav Hapla   PetscFunctionBegin;
5786a7c706bSVaclav Hapla   PetscSorted(n,X,*sorted);
5796a7c706bSVaclav Hapla   PetscFunctionReturn(0);
5806a7c706bSVaclav Hapla }
5816a7c706bSVaclav Hapla 
5826a7c706bSVaclav Hapla /*@
58317d7d925SStefano Zampini    PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order.
58417d7d925SStefano Zampini 
58517d7d925SStefano Zampini    Not Collective
58617d7d925SStefano Zampini 
58717d7d925SStefano Zampini    Input Parameters:
58817d7d925SStefano Zampini +  n  - number of values
58959e16d97SJunchao Zhang -  X  - array of integers
59017d7d925SStefano Zampini 
59117d7d925SStefano Zampini    Level: intermediate
59217d7d925SStefano Zampini 
593676f2a66SJacob Faibussowitsch    Notes:
594676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
595a5b23f4aSJose E. Roman    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
596676f2a66SJacob Faibussowitsch    code to see which routine is fastest.
597676f2a66SJacob Faibussowitsch 
598676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation()
59917d7d925SStefano Zampini @*/
600b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIInt(PetscInt n,PetscMPIInt X[])
60117d7d925SStefano Zampini {
60259e16d97SJunchao Zhang   PetscErrorCode ierr;
60359e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
60417d7d925SStefano Zampini 
60517d7d925SStefano Zampini   PetscFunctionBegin;
60659e16d97SJunchao Zhang   QuickSort1(PetscSortMPIInt,X,n,pivot,t1,ierr);
60717d7d925SStefano Zampini   PetscFunctionReturn(0);
60817d7d925SStefano Zampini }
60917d7d925SStefano Zampini 
61017d7d925SStefano Zampini /*@
61117d7d925SStefano Zampini    PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries
61217d7d925SStefano Zampini 
61317d7d925SStefano Zampini    Not Collective
61417d7d925SStefano Zampini 
61517d7d925SStefano Zampini    Input Parameters:
61617d7d925SStefano Zampini +  n  - number of values
61759e16d97SJunchao Zhang -  X  - array of integers
61817d7d925SStefano Zampini 
61917d7d925SStefano Zampini    Output Parameter:
62017d7d925SStefano Zampini .  n - number of non-redundant values
62117d7d925SStefano Zampini 
62217d7d925SStefano Zampini    Level: intermediate
62317d7d925SStefano Zampini 
62417d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
62517d7d925SStefano Zampini @*/
62659e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt X[])
62717d7d925SStefano Zampini {
62817d7d925SStefano Zampini   PetscErrorCode ierr;
62917d7d925SStefano Zampini   PetscInt       i,s = 0,N = *n, b = 0;
63017d7d925SStefano Zampini 
63117d7d925SStefano Zampini   PetscFunctionBegin;
63259e16d97SJunchao Zhang   ierr = PetscSortMPIInt(N,X);CHKERRQ(ierr);
63317d7d925SStefano Zampini   for (i=0; i<N-1; i++) {
63459e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
63559e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
636a297a907SKarl Rupp     } else s++;
63717d7d925SStefano Zampini   }
63817d7d925SStefano Zampini   *n = N - s;
63917d7d925SStefano Zampini   PetscFunctionReturn(0);
64017d7d925SStefano Zampini }
641c1f0200aSDmitry Karpeev 
6424d615ea0SBarry Smith /*@
6434d615ea0SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order;
6444d615ea0SBarry Smith        changes a second array to match the sorted first array.
6454d615ea0SBarry Smith 
6464d615ea0SBarry Smith    Not Collective
6474d615ea0SBarry Smith 
6484d615ea0SBarry Smith    Input Parameters:
6494d615ea0SBarry Smith +  n  - number of values
65059e16d97SJunchao Zhang .  X  - array of integers
65159e16d97SJunchao Zhang -  Y  - second array of integers
6524d615ea0SBarry Smith 
6534d615ea0SBarry Smith    Level: intermediate
6544d615ea0SBarry Smith 
65566f49146SPierre Jolivet .seealso: PetscMPIIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
6564d615ea0SBarry Smith @*/
657b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt X[],PetscMPIInt Y[])
6584d615ea0SBarry Smith {
6594d615ea0SBarry Smith   PetscErrorCode ierr;
66059e16d97SJunchao Zhang   PetscMPIInt    pivot,t1,t2;
6614d615ea0SBarry Smith 
6624d615ea0SBarry Smith   PetscFunctionBegin;
66359e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithArray,X,Y,n,pivot,t1,t2,ierr);
6645569a785SJunchao Zhang   PetscFunctionReturn(0);
6655569a785SJunchao Zhang }
6665569a785SJunchao Zhang 
6675569a785SJunchao Zhang /*@
6685569a785SJunchao Zhang    PetscSortMPIIntWithIntArray - Sorts an array of MPI integers in place in increasing order;
6695569a785SJunchao Zhang        changes a second array of Petsc intergers to match the sorted first array.
6705569a785SJunchao Zhang 
6715569a785SJunchao Zhang    Not Collective
6725569a785SJunchao Zhang 
6735569a785SJunchao Zhang    Input Parameters:
6745569a785SJunchao Zhang +  n  - number of values
67559e16d97SJunchao Zhang .  X  - array of MPI integers
67659e16d97SJunchao Zhang -  Y  - second array of Petsc integers
6775569a785SJunchao Zhang 
6785569a785SJunchao Zhang    Level: intermediate
6795569a785SJunchao Zhang 
6805569a785SJunchao Zhang    Notes: this routine is useful when one needs to sort MPI ranks with other integer arrays.
6815569a785SJunchao Zhang 
682676f2a66SJacob Faibussowitsch .seealso: PetscSortMPIIntWithArray(), PetscIntSortSemiOrderedWithArray(), PetscTimSortWithArray()
6835569a785SJunchao Zhang @*/
684b4f531b9SJunchao Zhang PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n,PetscMPIInt X[],PetscInt Y[])
6855569a785SJunchao Zhang {
6865569a785SJunchao Zhang   PetscErrorCode ierr;
68759e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
6885569a785SJunchao Zhang   PetscInt       t2;
6895569a785SJunchao Zhang 
6905569a785SJunchao Zhang   PetscFunctionBegin;
69159e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithIntArray,X,Y,n,pivot,t1,t2,ierr);
692e5c89e4eSSatish Balay   PetscFunctionReturn(0);
693e5c89e4eSSatish Balay }
694e5c89e4eSSatish Balay 
695e5c89e4eSSatish Balay /*@
696e5c89e4eSSatish Balay    PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
697e5c89e4eSSatish Balay        changes a second SCALAR array to match the sorted first INTEGER array.
698e5c89e4eSSatish Balay 
699e5c89e4eSSatish Balay    Not Collective
700e5c89e4eSSatish Balay 
701e5c89e4eSSatish Balay    Input Parameters:
702e5c89e4eSSatish Balay +  n  - number of values
70359e16d97SJunchao Zhang .  X  - array of integers
70459e16d97SJunchao Zhang -  Y  - second array of scalars
705e5c89e4eSSatish Balay 
706e5c89e4eSSatish Balay    Level: intermediate
707e5c89e4eSSatish Balay 
70866f49146SPierre Jolivet .seealso: PetscTimSortWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
709e5c89e4eSSatish Balay @*/
710b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithScalarArray(PetscInt n,PetscInt X[],PetscScalar Y[])
711e5c89e4eSSatish Balay {
712e5c89e4eSSatish Balay   PetscErrorCode ierr;
71359e16d97SJunchao Zhang   PetscInt       pivot,t1;
71459e16d97SJunchao Zhang   PetscScalar    t2;
715e5c89e4eSSatish Balay 
716e5c89e4eSSatish Balay   PetscFunctionBegin;
71759e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithScalarArray,X,Y,n,pivot,t1,t2,ierr);
71817df18f2SToby Isaac   PetscFunctionReturn(0);
71917df18f2SToby Isaac }
72017df18f2SToby Isaac 
7216c2863d0SBarry Smith /*@C
72217df18f2SToby Isaac    PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order;
72317df18f2SToby Isaac        changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
72417df18f2SToby Isaac        provide workspace (the size of an element in the data array) to use when sorting.
72517df18f2SToby Isaac 
72617df18f2SToby Isaac    Not Collective
72717df18f2SToby Isaac 
72817df18f2SToby Isaac    Input Parameters:
72917df18f2SToby Isaac +  n  - number of values
73059e16d97SJunchao Zhang .  X  - array of integers
73159e16d97SJunchao Zhang .  Y  - second array of data
73217df18f2SToby Isaac .  size - sizeof elements in the data array in bytes
73359e16d97SJunchao Zhang -  t2   - workspace of "size" bytes used when sorting
73417df18f2SToby Isaac 
73517df18f2SToby Isaac    Level: intermediate
73617df18f2SToby Isaac 
73766f49146SPierre Jolivet .seealso: PetscTimSortWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
73817df18f2SToby Isaac @*/
739b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithDataArray(PetscInt n,PetscInt X[],void *Y,size_t size,void *t2)
74017df18f2SToby Isaac {
74117df18f2SToby Isaac   PetscErrorCode ierr;
74259e16d97SJunchao Zhang   char           *YY = (char*)Y;
74359e16d97SJunchao Zhang   PetscInt       i,j,p,t1,pivot,hi=n-1,l,r;
74417df18f2SToby Isaac 
74517df18f2SToby Isaac   PetscFunctionBegin;
74617df18f2SToby Isaac   if (n<8) {
74759e16d97SJunchao Zhang     for (i=0; i<n; i++) {
74859e16d97SJunchao Zhang       pivot = X[i];
74959e16d97SJunchao Zhang       for (j=i+1; j<n; j++) {
75059e16d97SJunchao Zhang         if (pivot > X[j]) {
75159e16d97SJunchao Zhang           SWAP2Data(X[i],X[j],YY+size*i,YY+size*j,t1,t2,size);
75259e16d97SJunchao Zhang           pivot = X[i];
75317df18f2SToby Isaac         }
75417df18f2SToby Isaac       }
75517df18f2SToby Isaac     }
75617df18f2SToby Isaac   } else {
75759e16d97SJunchao Zhang     /* Two way partition */
75859e16d97SJunchao Zhang     p     = MEDIAN(X,hi);
75959e16d97SJunchao Zhang     pivot = X[p];
76059e16d97SJunchao Zhang     l     = 0;
76159e16d97SJunchao Zhang     r     = hi;
76259e16d97SJunchao Zhang     while (1) {
76359e16d97SJunchao Zhang       while (X[l] < pivot) l++;
76459e16d97SJunchao Zhang       while (X[r] > pivot) r--;
76559e16d97SJunchao Zhang       if (l >= r) {r++; break;}
76659e16d97SJunchao Zhang       SWAP2Data(X[l],X[r],YY+size*l,YY+size*r,t1,t2,size);
76759e16d97SJunchao Zhang       l++;
76859e16d97SJunchao Zhang       r--;
76959e16d97SJunchao Zhang     }
77059e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(l,X,Y,size,t2);CHKERRQ(ierr);
77159e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(hi-r+1,X+r,YY+size*r,size,t2);CHKERRQ(ierr);
77217df18f2SToby Isaac   }
77317df18f2SToby Isaac   PetscFunctionReturn(0);
77417df18f2SToby Isaac }
77517df18f2SToby Isaac 
77621e72a00SBarry Smith /*@
77721e72a00SBarry Smith    PetscMergeIntArray -     Merges two SORTED integer arrays, removes duplicate elements.
77821e72a00SBarry Smith 
77921e72a00SBarry Smith    Not Collective
78021e72a00SBarry Smith 
78121e72a00SBarry Smith    Input Parameters:
78221e72a00SBarry Smith +  an  - number of values in the first array
78321e72a00SBarry Smith .  aI  - first sorted array of integers
78421e72a00SBarry Smith .  bn  - number of values in the second array
78521e72a00SBarry Smith -  bI  - second array of integers
78621e72a00SBarry Smith 
78721e72a00SBarry Smith    Output Parameters:
78821e72a00SBarry Smith +  n   - number of values in the merged array
7896c2863d0SBarry Smith -  L   - merged sorted array, this is allocated if an array is not provided
79021e72a00SBarry Smith 
79121e72a00SBarry Smith    Level: intermediate
79221e72a00SBarry Smith 
79366f49146SPierre Jolivet .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
79421e72a00SBarry Smith @*/
7956c2863d0SBarry Smith PetscErrorCode  PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
79621e72a00SBarry Smith {
79721e72a00SBarry Smith   PetscErrorCode ierr;
79821e72a00SBarry Smith   PetscInt       *L_ = *L, ak, bk, k;
79921e72a00SBarry Smith 
800362febeeSStefano Zampini   PetscFunctionBegin;
80121e72a00SBarry Smith   if (!L_) {
80221e72a00SBarry Smith     ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr);
80321e72a00SBarry Smith     L_   = *L;
80421e72a00SBarry Smith   }
80521e72a00SBarry Smith   k = ak = bk = 0;
80621e72a00SBarry Smith   while (ak < an && bk < bn) {
80721e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
80821e72a00SBarry Smith       L_[k] = aI[ak];
80921e72a00SBarry Smith       ++ak;
81021e72a00SBarry Smith       ++bk;
81121e72a00SBarry Smith       ++k;
81221e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
81321e72a00SBarry Smith       L_[k] = aI[ak];
81421e72a00SBarry Smith       ++ak;
81521e72a00SBarry Smith       ++k;
81621e72a00SBarry Smith     } else {
81721e72a00SBarry Smith       L_[k] = bI[bk];
81821e72a00SBarry Smith       ++bk;
81921e72a00SBarry Smith       ++k;
82021e72a00SBarry Smith     }
82121e72a00SBarry Smith   }
82221e72a00SBarry Smith   if (ak < an) {
823580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
82421e72a00SBarry Smith     k   += (an-ak);
82521e72a00SBarry Smith   }
82621e72a00SBarry Smith   if (bk < bn) {
827580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
82821e72a00SBarry Smith     k   += (bn-bk);
82921e72a00SBarry Smith   }
83021e72a00SBarry Smith   *n = k;
83121e72a00SBarry Smith   PetscFunctionReturn(0);
83221e72a00SBarry Smith }
833b4301105SBarry Smith 
834c1f0200aSDmitry Karpeev /*@
83521e72a00SBarry Smith    PetscMergeIntArrayPair -     Merges two SORTED integer arrays that share NO common values along with an additional array of integers.
836c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
837c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
838c1f0200aSDmitry Karpeev 
839c1f0200aSDmitry Karpeev    Not Collective
840c1f0200aSDmitry Karpeev 
841c1f0200aSDmitry Karpeev    Input Parameters:
842c1f0200aSDmitry Karpeev +  an  - number of values in the first array
843c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
844c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
845c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
846c1f0200aSDmitry Karpeev .  bI  - second array of integers
847c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
848c1f0200aSDmitry Karpeev 
849c1f0200aSDmitry Karpeev    Output Parameters:
850c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
85114c49006SJed Brown .  L   - merged sorted array
852c1f0200aSDmitry Karpeev -  J   - merged additional array
853c1f0200aSDmitry Karpeev 
85495452b02SPatrick Sanan    Notes:
85595452b02SPatrick 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
856c1f0200aSDmitry Karpeev    Level: intermediate
857c1f0200aSDmitry Karpeev 
85866f49146SPierre Jolivet .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
859c1f0200aSDmitry Karpeev @*/
8606c2863d0SBarry 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)
861c1f0200aSDmitry Karpeev {
862c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
86370d8d27cSBarry Smith   PetscInt       n_, *L_, *J_, ak, bk, k;
864c1f0200aSDmitry Karpeev 
86514c49006SJed Brown   PetscFunctionBegin;
86670d8d27cSBarry Smith   PetscValidIntPointer(L,8);
86770d8d27cSBarry Smith   PetscValidIntPointer(J,9);
868c1f0200aSDmitry Karpeev   n_ = an + bn;
869c1f0200aSDmitry Karpeev   *n = n_;
87070d8d27cSBarry Smith   if (!*L) {
871785e854fSJed Brown     ierr = PetscMalloc1(n_, L);CHKERRQ(ierr);
87270d8d27cSBarry Smith   }
873d7aa01a8SHong Zhang   L_ = *L;
87470d8d27cSBarry Smith   if (!*J) {
87570d8d27cSBarry Smith     ierr = PetscMalloc1(n_, J);CHKERRQ(ierr);
876c1f0200aSDmitry Karpeev   }
877c1f0200aSDmitry Karpeev   J_   = *J;
878c1f0200aSDmitry Karpeev   k = ak = bk = 0;
879c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
880c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
881d7aa01a8SHong Zhang       L_[k] = aI[ak];
882c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
883c1f0200aSDmitry Karpeev       ++ak;
884c1f0200aSDmitry Karpeev       ++k;
885a297a907SKarl Rupp     } else {
886d7aa01a8SHong Zhang       L_[k] = bI[bk];
887c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
888c1f0200aSDmitry Karpeev       ++bk;
889c1f0200aSDmitry Karpeev       ++k;
890c1f0200aSDmitry Karpeev     }
891c1f0200aSDmitry Karpeev   }
892c1f0200aSDmitry Karpeev   if (ak < an) {
893580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
894580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,aJ+ak,an-ak);CHKERRQ(ierr);
895c1f0200aSDmitry Karpeev     k   += (an-ak);
896c1f0200aSDmitry Karpeev   }
897c1f0200aSDmitry Karpeev   if (bk < bn) {
898580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
899580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,bJ+bk,bn-bk);CHKERRQ(ierr);
900c1f0200aSDmitry Karpeev   }
901c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
902c1f0200aSDmitry Karpeev }
903c1f0200aSDmitry Karpeev 
904e498c390SJed Brown /*@
905e498c390SJed Brown    PetscMergeMPIIntArray -     Merges two SORTED integer arrays.
906e498c390SJed Brown 
907e498c390SJed Brown    Not Collective
908e498c390SJed Brown 
909e498c390SJed Brown    Input Parameters:
910e498c390SJed Brown +  an  - number of values in the first array
911e498c390SJed Brown .  aI  - first sorted array of integers
912e498c390SJed Brown .  bn  - number of values in the second array
913e498c390SJed Brown -  bI  - second array of integers
914e498c390SJed Brown 
915e498c390SJed Brown    Output Parameters:
916e498c390SJed Brown +  n   - number of values in the merged array (<= an + bn)
917e498c390SJed Brown -  L   - merged sorted array, allocated if address of NULL pointer is passed
918e498c390SJed Brown 
919e498c390SJed Brown    Level: intermediate
920e498c390SJed Brown 
92166f49146SPierre Jolivet .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
922e498c390SJed Brown @*/
923e498c390SJed Brown PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L)
924e498c390SJed Brown {
925e498c390SJed Brown   PetscErrorCode ierr;
926e498c390SJed Brown   PetscInt       ai,bi,k;
927e498c390SJed Brown 
928e498c390SJed Brown   PetscFunctionBegin;
929e498c390SJed Brown   if (!*L) {ierr = PetscMalloc1((an+bn),L);CHKERRQ(ierr);}
930e498c390SJed Brown   for (ai=0,bi=0,k=0; ai<an || bi<bn;) {
931e498c390SJed Brown     PetscInt t = -1;
932e498c390SJed Brown     for (; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
933e498c390SJed Brown     for (; bi<bn && bI[bi] == t; bi++);
934e498c390SJed Brown     for (; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
935e498c390SJed Brown     for (; ai<an && aI[ai] == t; ai++);
936e498c390SJed Brown   }
937e498c390SJed Brown   *n = k;
938e498c390SJed Brown   PetscFunctionReturn(0);
939e498c390SJed Brown }
940e5c89e4eSSatish Balay 
9416c2863d0SBarry Smith /*@C
94242eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
94342eaa7f3SBarry Smith 
94442eaa7f3SBarry Smith    Not Collective
94542eaa7f3SBarry Smith 
94642eaa7f3SBarry Smith    Input Parameters:
94742eaa7f3SBarry Smith +  n  - number of values
94842eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
94942eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
95042eaa7f3SBarry Smith 
95142eaa7f3SBarry Smith    Output Parameters:
95242eaa7f3SBarry Smith +  Nlevels - the number of levels
95342eaa7f3SBarry Smith .  Level - for each node tells its level
95442eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
95542eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
95642eaa7f3SBarry Smith -  Column - for each id tells its column index
95742eaa7f3SBarry Smith 
9586c2863d0SBarry Smith    Level: developer
95942eaa7f3SBarry Smith 
96095452b02SPatrick Sanan    Notes:
96195452b02SPatrick Sanan     This code is not currently used
96242eaa7f3SBarry Smith 
96342eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation()
96442eaa7f3SBarry Smith @*/
9657087cfbeSBarry Smith PetscErrorCode  PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column)
96642eaa7f3SBarry Smith {
96742eaa7f3SBarry Smith   PetscInt       i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column;
96842eaa7f3SBarry Smith   PetscErrorCode ierr;
969ace3abfcSBarry Smith   PetscBool      done = PETSC_FALSE;
97042eaa7f3SBarry Smith 
97142eaa7f3SBarry Smith   PetscFunctionBegin;
9722c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mask[0],PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set");
97342eaa7f3SBarry Smith   for (i=0; i<n; i++) {
97442eaa7f3SBarry Smith     if (mask[i]) continue;
9752c71b3e2SJacob Faibussowitsch     PetscCheckFalse(parentid[i]  == i,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent");
9762c71b3e2SJacob Faibussowitsch     PetscCheckFalse(parentid[i] && mask[parentid[i]],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked");
97742eaa7f3SBarry Smith   }
97842eaa7f3SBarry Smith 
97942eaa7f3SBarry Smith   for (i=0; i<n; i++) {
98042eaa7f3SBarry Smith     if (!mask[i]) nmask++;
98142eaa7f3SBarry Smith   }
98242eaa7f3SBarry Smith 
98342eaa7f3SBarry Smith   /* determine the level in the tree of each node */
9841795a4d1SJed Brown   ierr = PetscCalloc1(n,&level);CHKERRQ(ierr);
985a297a907SKarl Rupp 
98642eaa7f3SBarry Smith   level[0] = 1;
98742eaa7f3SBarry Smith   while (!done) {
98842eaa7f3SBarry Smith     done = PETSC_TRUE;
98942eaa7f3SBarry Smith     for (i=0; i<n; i++) {
99042eaa7f3SBarry Smith       if (mask[i]) continue;
99142eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
99242eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
99342eaa7f3SBarry Smith     }
99442eaa7f3SBarry Smith   }
99542eaa7f3SBarry Smith   for (i=0; i<n; i++) {
99642eaa7f3SBarry Smith     level[i]--;
99742eaa7f3SBarry Smith     nlevels = PetscMax(nlevels,level[i]);
99842eaa7f3SBarry Smith   }
99942eaa7f3SBarry Smith 
100042eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
10011795a4d1SJed Brown   ierr = PetscCalloc1(nlevels,&levelcnt);CHKERRQ(ierr);
100242eaa7f3SBarry Smith   for (i=0; i<n; i++) {
100342eaa7f3SBarry Smith     if (mask[i]) continue;
100442eaa7f3SBarry Smith     levelcnt[level[i]-1]++;
100542eaa7f3SBarry Smith   }
1006a297a907SKarl Rupp   for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]);
100742eaa7f3SBarry Smith 
100842eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
1009dcca6d9dSJed Brown   ierr = PetscMalloc2(levelmax,&workid,levelmax,&workparentid);CHKERRQ(ierr);
1010785e854fSJed Brown   ierr = PetscMalloc1(nmask,&idbylevel);CHKERRQ(ierr);
101142eaa7f3SBarry Smith   for (j=1; j<=nlevels;j++) {
101242eaa7f3SBarry Smith     cnt = 0;
101342eaa7f3SBarry Smith     for (i=0; i<n; i++) {
101442eaa7f3SBarry Smith       if (mask[i]) continue;
101542eaa7f3SBarry Smith       if (level[i] != j) continue;
101642eaa7f3SBarry Smith       workid[cnt]         = i;
101742eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
101842eaa7f3SBarry Smith     }
101942eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
102042eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
102142eaa7f3SBarry Smith     ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr);
102242eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
102342eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
1024580bdb30SBarry Smith     ierr  = PetscArraycpy(idbylevel+tcnt,workid,cnt);CHKERRQ(ierr);
102542eaa7f3SBarry Smith     tcnt += cnt;
102642eaa7f3SBarry Smith   }
10272c71b3e2SJacob Faibussowitsch   PetscCheckFalse(tcnt != nmask,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes");
102842eaa7f3SBarry Smith   ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr);
102942eaa7f3SBarry Smith 
103042eaa7f3SBarry Smith   /* for each node list its column */
1031785e854fSJed Brown   ierr = PetscMalloc1(n,&column);CHKERRQ(ierr);
103242eaa7f3SBarry Smith   cnt = 0;
103342eaa7f3SBarry Smith   for (j=0; j<nlevels; j++) {
103442eaa7f3SBarry Smith     for (i=0; i<levelcnt[j]; i++) {
103542eaa7f3SBarry Smith       column[idbylevel[cnt++]] = i;
103642eaa7f3SBarry Smith     }
103742eaa7f3SBarry Smith   }
103842eaa7f3SBarry Smith 
103942eaa7f3SBarry Smith   *Nlevels   = nlevels;
104042eaa7f3SBarry Smith   *Level     = level;
104142eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
104242eaa7f3SBarry Smith   *Idbylevel = idbylevel;
104342eaa7f3SBarry Smith   *Column    = column;
104442eaa7f3SBarry Smith   PetscFunctionReturn(0);
104542eaa7f3SBarry Smith }
1046ce605777SToby Isaac 
1047ce605777SToby Isaac /*@
1048ce605777SToby Isaac   PetscParallelSortedInt - Check whether an integer array, distributed over a communicator, is globally sorted.
1049ce605777SToby Isaac 
1050ce605777SToby Isaac   Collective
1051ce605777SToby Isaac 
1052ce605777SToby Isaac   Input Parameters:
1053ce605777SToby Isaac + comm - the MPI communicator
1054ce605777SToby Isaac . n - the local number of integers
1055ce605777SToby Isaac - keys - the local array of integers
1056ce605777SToby Isaac 
1057ce605777SToby Isaac   Output Parameters:
1058ce605777SToby Isaac . is_sorted - whether the array is globally sorted
1059ce605777SToby Isaac 
1060ce605777SToby Isaac   Level: developer
1061ce605777SToby Isaac 
1062ce605777SToby Isaac .seealso: PetscParallelSortInt()
1063ce605777SToby Isaac @*/
1064ce605777SToby Isaac PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
1065ce605777SToby Isaac {
1066ce605777SToby Isaac   PetscBool      sorted;
1067ce605777SToby Isaac   PetscInt       i, min, max, prevmax;
10682a1da528SToby Isaac   PetscMPIInt    rank;
1069ce605777SToby Isaac   PetscErrorCode ierr;
1070ce605777SToby Isaac 
1071ce605777SToby Isaac   PetscFunctionBegin;
1072ce605777SToby Isaac   sorted = PETSC_TRUE;
1073ce605777SToby Isaac   min    = PETSC_MAX_INT;
1074ce605777SToby Isaac   max    = PETSC_MIN_INT;
1075ce605777SToby Isaac   if (n) {
1076ce605777SToby Isaac     min = keys[0];
1077ce605777SToby Isaac     max = keys[0];
1078ce605777SToby Isaac   }
1079ce605777SToby Isaac   for (i = 1; i < n; i++) {
1080ce605777SToby Isaac     if (keys[i] < keys[i - 1]) break;
1081ce605777SToby Isaac     min = PetscMin(min,keys[i]);
1082ce605777SToby Isaac     max = PetscMax(max,keys[i]);
1083ce605777SToby Isaac   }
1084ce605777SToby Isaac   if (i < n) sorted = PETSC_FALSE;
1085ce605777SToby Isaac   prevmax = PETSC_MIN_INT;
1086ffc4695bSBarry Smith   ierr = MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(ierr);
1087ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
1088dd400576SPatrick Sanan   if (rank == 0) prevmax = PETSC_MIN_INT;
1089ce605777SToby Isaac   if (prevmax > min) sorted = PETSC_FALSE;
1090ffc4695bSBarry Smith   ierr = MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRMPI(ierr);
1091ce605777SToby Isaac   PetscFunctionReturn(0);
1092ce605777SToby Isaac }
1093