xref: /petsc/src/sys/utils/sorti.c (revision b6a18d215e416d140752c5e92bb1511956dfe2a5)
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 {                                                                           \
126981bb840SJunchao 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 {                                                                           \
149981bb840SJunchao 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 {                                                                           \
171981bb840SJunchao 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 {                                                                           \
193981bb840SJunchao 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 /*@
319*b6a18d21SVaclav Hapla    PetscSortedCheckDupsInt - Checks if a sorted integer array has duplicates
320*b6a18d21SVaclav Hapla 
321*b6a18d21SVaclav Hapla    Not Collective
322*b6a18d21SVaclav Hapla 
323*b6a18d21SVaclav Hapla    Input Parameters:
324*b6a18d21SVaclav Hapla +  n  - number of values
325*b6a18d21SVaclav Hapla -  X  - sorted array of integers
326*b6a18d21SVaclav Hapla 
327*b6a18d21SVaclav Hapla    Output Parameter:
328*b6a18d21SVaclav Hapla .  dups - True if the array has dups, otherwise false
329*b6a18d21SVaclav Hapla 
330*b6a18d21SVaclav Hapla    Level: intermediate
331*b6a18d21SVaclav Hapla 
332*b6a18d21SVaclav Hapla .seealso: PetscSortInt(), PetscCheckDupsInt(), PetscSortRemoveDupsInt(), PetscSortedRemoveDupsInt()
333*b6a18d21SVaclav Hapla @*/
334*b6a18d21SVaclav Hapla PetscErrorCode  PetscSortedCheckDupsInt(PetscInt n,const PetscInt X[],PetscBool *flg)
335*b6a18d21SVaclav Hapla {
336*b6a18d21SVaclav Hapla   PetscInt i;
337*b6a18d21SVaclav Hapla 
338*b6a18d21SVaclav Hapla   PetscFunctionBegin;
339*b6a18d21SVaclav Hapla   PetscCheckSorted(n,X);
340*b6a18d21SVaclav Hapla   *flg = PETSC_FALSE;
341*b6a18d21SVaclav Hapla   for (i=0; i<n-1; i++) {
342*b6a18d21SVaclav Hapla     if (X[i+1] == X[i]) {
343*b6a18d21SVaclav Hapla       *flg = PETSC_TRUE;
344*b6a18d21SVaclav Hapla       break;
345*b6a18d21SVaclav Hapla     }
346*b6a18d21SVaclav Hapla   }
347*b6a18d21SVaclav Hapla   PetscFunctionReturn(0);
348*b6a18d21SVaclav Hapla }
349*b6a18d21SVaclav Hapla 
350*b6a18d21SVaclav Hapla /*@
3511db36a52SBarry Smith    PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries
3521db36a52SBarry Smith 
3531db36a52SBarry Smith    Not Collective
3541db36a52SBarry Smith 
3551db36a52SBarry Smith    Input Parameters:
3561db36a52SBarry Smith +  n  - number of values
35759e16d97SJunchao Zhang -  X  - array of integers
3581db36a52SBarry Smith 
3591db36a52SBarry Smith    Output Parameter:
3601db36a52SBarry Smith .  n - number of non-redundant values
3611db36a52SBarry Smith 
3621db36a52SBarry Smith    Level: intermediate
3631db36a52SBarry Smith 
364676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortedRemoveDupsInt()
3651db36a52SBarry Smith @*/
36659e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsInt(PetscInt *n,PetscInt X[])
3671db36a52SBarry Smith {
3681db36a52SBarry Smith   PetscErrorCode ierr;
3691db36a52SBarry Smith 
3701db36a52SBarry Smith   PetscFunctionBegin;
37159e16d97SJunchao Zhang   ierr = PetscSortInt(*n,X);CHKERRQ(ierr);
37259e16d97SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(n,X);CHKERRQ(ierr);
3731db36a52SBarry Smith   PetscFunctionReturn(0);
3741db36a52SBarry Smith }
3751db36a52SBarry Smith 
37660e03357SMatthew G Knepley /*@
377d9f1162dSJed Brown   PetscFindInt - Finds integer in a sorted array of integers
37860e03357SMatthew G Knepley 
37960e03357SMatthew G Knepley    Not Collective
38060e03357SMatthew G Knepley 
38160e03357SMatthew G Knepley    Input Parameters:
38260e03357SMatthew G Knepley +  key - the integer to locate
383d9f1162dSJed Brown .  n   - number of values in the array
38459e16d97SJunchao Zhang -  X  - array of integers
38560e03357SMatthew G Knepley 
38660e03357SMatthew G Knepley    Output Parameter:
387d9f1162dSJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
38860e03357SMatthew G Knepley 
38960e03357SMatthew G Knepley    Level: intermediate
39060e03357SMatthew G Knepley 
391676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
39260e03357SMatthew G Knepley @*/
39359e16d97SJunchao Zhang PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
39460e03357SMatthew G Knepley {
395d9f1162dSJed Brown   PetscInt lo = 0,hi = n;
39660e03357SMatthew G Knepley 
39760e03357SMatthew G Knepley   PetscFunctionBegin;
39860e03357SMatthew G Knepley   PetscValidPointer(loc,4);
39998781d82SMatthew G Knepley   if (!n) {*loc = -1; PetscFunctionReturn(0);}
40059e16d97SJunchao Zhang   PetscValidPointer(X,3);
4016a7c706bSVaclav Hapla   PetscCheckSorted(n,X);
402d9f1162dSJed Brown   while (hi - lo > 1) {
403d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo)/2;
40459e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
405d9f1162dSJed Brown     else               lo = mid;
40660e03357SMatthew G Knepley   }
40759e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
40860e03357SMatthew G Knepley   PetscFunctionReturn(0);
40960e03357SMatthew G Knepley }
41060e03357SMatthew G Knepley 
411d2aeb606SJed Brown /*@
412f1cab4e1SJunchao Zhang   PetscCheckDupsInt - Checks if an integer array has duplicates
413f1cab4e1SJunchao Zhang 
414f1cab4e1SJunchao Zhang    Not Collective
415f1cab4e1SJunchao Zhang 
416f1cab4e1SJunchao Zhang    Input Parameters:
417f1cab4e1SJunchao Zhang +  n  - number of values in the array
418f1cab4e1SJunchao Zhang -  X  - array of integers
419f1cab4e1SJunchao Zhang 
420f1cab4e1SJunchao Zhang    Output Parameter:
421f1cab4e1SJunchao Zhang .  dups - True if the array has dups, otherwise false
422f1cab4e1SJunchao Zhang 
423f1cab4e1SJunchao Zhang    Level: intermediate
424f1cab4e1SJunchao Zhang 
425*b6a18d21SVaclav Hapla .seealso: PetscSortRemoveDupsInt(), PetscSortedCheckDupsInt()
426f1cab4e1SJunchao Zhang @*/
427f1cab4e1SJunchao Zhang PetscErrorCode PetscCheckDupsInt(PetscInt n,const PetscInt X[],PetscBool *dups)
428f1cab4e1SJunchao Zhang {
429f1cab4e1SJunchao Zhang   PetscErrorCode ierr;
430f1cab4e1SJunchao Zhang   PetscInt       i;
431f1cab4e1SJunchao Zhang   PetscHSetI     ht;
432f1cab4e1SJunchao Zhang   PetscBool      missing;
433f1cab4e1SJunchao Zhang 
434f1cab4e1SJunchao Zhang   PetscFunctionBegin;
435f1cab4e1SJunchao Zhang   PetscValidPointer(dups,3);
436f1cab4e1SJunchao Zhang   *dups = PETSC_FALSE;
437f1cab4e1SJunchao Zhang   if (n > 1) {
438f1cab4e1SJunchao Zhang     ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
439f1cab4e1SJunchao Zhang     ierr = PetscHSetIResize(ht,n);CHKERRQ(ierr);
440f1cab4e1SJunchao Zhang     for (i=0; i<n; i++) {
441f1cab4e1SJunchao Zhang       ierr = PetscHSetIQueryAdd(ht,X[i],&missing);CHKERRQ(ierr);
442f1cab4e1SJunchao Zhang       if (!missing) {*dups = PETSC_TRUE; break;}
443f1cab4e1SJunchao Zhang     }
444f1cab4e1SJunchao Zhang     ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
445f1cab4e1SJunchao Zhang   }
446f1cab4e1SJunchao Zhang   PetscFunctionReturn(0);
447f1cab4e1SJunchao Zhang }
448f1cab4e1SJunchao Zhang 
449f1cab4e1SJunchao Zhang /*@
450d2aeb606SJed Brown   PetscFindMPIInt - Finds MPI integer in a sorted array of integers
451d2aeb606SJed Brown 
452d2aeb606SJed Brown    Not Collective
453d2aeb606SJed Brown 
454d2aeb606SJed Brown    Input Parameters:
455d2aeb606SJed Brown +  key - the integer to locate
456d2aeb606SJed Brown .  n   - number of values in the array
45759e16d97SJunchao Zhang -  X   - array of integers
458d2aeb606SJed Brown 
459d2aeb606SJed Brown    Output Parameter:
460d2aeb606SJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
461d2aeb606SJed Brown 
462d2aeb606SJed Brown    Level: intermediate
463d2aeb606SJed Brown 
464676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
465d2aeb606SJed Brown @*/
46659e16d97SJunchao Zhang PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
467d2aeb606SJed Brown {
468d2aeb606SJed Brown   PetscInt lo = 0,hi = n;
469d2aeb606SJed Brown 
470d2aeb606SJed Brown   PetscFunctionBegin;
471d2aeb606SJed Brown   PetscValidPointer(loc,4);
472d2aeb606SJed Brown   if (!n) {*loc = -1; PetscFunctionReturn(0);}
47359e16d97SJunchao Zhang   PetscValidPointer(X,3);
4746a7c706bSVaclav Hapla   PetscCheckSorted(n,X);
475d2aeb606SJed Brown   while (hi - lo > 1) {
476d2aeb606SJed Brown     PetscInt mid = lo + (hi - lo)/2;
47759e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
478d2aeb606SJed Brown     else               lo = mid;
479d2aeb606SJed Brown   }
48059e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
481e5c89e4eSSatish Balay   PetscFunctionReturn(0);
482e5c89e4eSSatish Balay }
483e5c89e4eSSatish Balay 
484e5c89e4eSSatish Balay /*@
485e5c89e4eSSatish Balay    PetscSortIntWithArray - Sorts an array of integers in place in increasing order;
486e5c89e4eSSatish Balay        changes a second array to match the sorted first array.
487e5c89e4eSSatish Balay 
488e5c89e4eSSatish Balay    Not Collective
489e5c89e4eSSatish Balay 
490e5c89e4eSSatish Balay    Input Parameters:
491e5c89e4eSSatish Balay +  n  - number of values
49259e16d97SJunchao Zhang .  X  - array of integers
49359e16d97SJunchao Zhang -  Y  - second array of integers
494e5c89e4eSSatish Balay 
495e5c89e4eSSatish Balay    Level: intermediate
496e5c89e4eSSatish Balay 
497981bb840SJunchao Zhang .seealso: PetscIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithCountArray()
498e5c89e4eSSatish Balay @*/
499b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithArray(PetscInt n,PetscInt X[],PetscInt Y[])
500e5c89e4eSSatish Balay {
501e5c89e4eSSatish Balay   PetscErrorCode ierr;
50259e16d97SJunchao Zhang   PetscInt       pivot,t1,t2;
503e5c89e4eSSatish Balay 
504e5c89e4eSSatish Balay   PetscFunctionBegin;
50559e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithArray,X,Y,n,pivot,t1,t2,ierr);
506c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
507c1f0200aSDmitry Karpeev }
508c1f0200aSDmitry Karpeev 
509c1f0200aSDmitry Karpeev /*@
510c1f0200aSDmitry Karpeev    PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order;
511c1f0200aSDmitry Karpeev        changes a pair of integer arrays to match the sorted first array.
512c1f0200aSDmitry Karpeev 
513c1f0200aSDmitry Karpeev    Not Collective
514c1f0200aSDmitry Karpeev 
515c1f0200aSDmitry Karpeev    Input Parameters:
516c1f0200aSDmitry Karpeev +  n  - number of values
51759e16d97SJunchao Zhang .  X  - array of integers
51859e16d97SJunchao Zhang .  Y  - second array of integers (first array of the pair)
51959e16d97SJunchao Zhang -  Z  - third array of integers  (second array of the pair)
520c1f0200aSDmitry Karpeev 
521c1f0200aSDmitry Karpeev    Level: intermediate
522c1f0200aSDmitry Karpeev 
523981bb840SJunchao Zhang .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortIntWithArray(), PetscIntSortSemiOrdered(), PetscSortIntWithIntCountArrayPair()
524c1f0200aSDmitry Karpeev @*/
525b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithArrayPair(PetscInt n,PetscInt X[],PetscInt Y[],PetscInt Z[])
526c1f0200aSDmitry Karpeev {
527c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
52859e16d97SJunchao Zhang   PetscInt       pivot,t1,t2,t3;
529c1f0200aSDmitry Karpeev 
530c1f0200aSDmitry Karpeev   PetscFunctionBegin;
53159e16d97SJunchao Zhang   QuickSort3(PetscSortIntWithArrayPair,X,Y,Z,n,pivot,t1,t2,t3,ierr);
532c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
533c1f0200aSDmitry Karpeev }
534c1f0200aSDmitry Karpeev 
53517d7d925SStefano Zampini /*@
536981bb840SJunchao Zhang    PetscSortIntWithCountArray - Sorts an array of integers in place in increasing order;
537981bb840SJunchao Zhang        changes a second array to match the sorted first array.
538981bb840SJunchao Zhang 
539981bb840SJunchao Zhang    Not Collective
540981bb840SJunchao Zhang 
541981bb840SJunchao Zhang    Input Parameters:
542981bb840SJunchao Zhang +  n  - number of values
543981bb840SJunchao Zhang .  X  - array of integers
544981bb840SJunchao Zhang -  Y  - second array of PetscCounts (signed integers)
545981bb840SJunchao Zhang 
546981bb840SJunchao Zhang    Level: intermediate
547981bb840SJunchao Zhang 
548981bb840SJunchao Zhang .seealso: PetscIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
549981bb840SJunchao Zhang @*/
550981bb840SJunchao Zhang PetscErrorCode  PetscSortIntWithCountArray(PetscCount n,PetscInt X[],PetscCount Y[])
551981bb840SJunchao Zhang {
552981bb840SJunchao Zhang   PetscErrorCode ierr;
553981bb840SJunchao Zhang   PetscInt       pivot,t1;
554981bb840SJunchao Zhang   PetscCount     t2;
555981bb840SJunchao Zhang 
556981bb840SJunchao Zhang   PetscFunctionBegin;
557981bb840SJunchao Zhang   QuickSort2(PetscSortIntWithCountArray,X,Y,n,pivot,t1,t2,ierr);
558981bb840SJunchao Zhang   PetscFunctionReturn(0);
559981bb840SJunchao Zhang }
560981bb840SJunchao Zhang 
561981bb840SJunchao Zhang /*@
562981bb840SJunchao Zhang    PetscSortIntWithIntCountArrayPair - Sorts an array of integers in place in increasing order;
563981bb840SJunchao Zhang        changes an integer array and a PetscCount array to match the sorted first array.
564981bb840SJunchao Zhang 
565981bb840SJunchao Zhang    Not Collective
566981bb840SJunchao Zhang 
567981bb840SJunchao Zhang    Input Parameters:
568981bb840SJunchao Zhang +  n  - number of values
569981bb840SJunchao Zhang .  X  - array of integers
570981bb840SJunchao Zhang .  Y  - second array of integers (first array of the pair)
571981bb840SJunchao Zhang -  Z  - third array of PetscCounts  (second array of the pair)
572981bb840SJunchao Zhang 
573981bb840SJunchao Zhang    Level: intermediate
574981bb840SJunchao Zhang 
575981bb840SJunchao Zhang    Notes:
576981bb840SJunchao 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.
577981bb840SJunchao Zhang 
578981bb840SJunchao Zhang .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray(), PetscIntSortSemiOrdered(), PetscSortIntWithArrayPair()
579981bb840SJunchao Zhang @*/
580981bb840SJunchao Zhang PetscErrorCode  PetscSortIntWithIntCountArrayPair(PetscCount n,PetscInt X[],PetscInt Y[],PetscCount Z[])
581981bb840SJunchao Zhang {
582981bb840SJunchao Zhang   PetscErrorCode ierr;
583981bb840SJunchao Zhang   PetscInt       pivot,t1,t2; /* pivot is take from X[], so its type is still PetscInt */
584981bb840SJunchao Zhang   PetscCount     t3; /* temp for Z[] */
585981bb840SJunchao Zhang 
586981bb840SJunchao Zhang   PetscFunctionBegin;
587981bb840SJunchao Zhang   QuickSort3(PetscSortIntWithIntCountArrayPair,X,Y,Z,n,pivot,t1,t2,t3,ierr);
588981bb840SJunchao Zhang   PetscFunctionReturn(0);
589981bb840SJunchao Zhang }
590981bb840SJunchao Zhang 
591981bb840SJunchao Zhang /*@
5926a7c706bSVaclav Hapla    PetscSortedMPIInt - Determines whether the array is sorted.
5936a7c706bSVaclav Hapla 
5946a7c706bSVaclav Hapla    Not Collective
5956a7c706bSVaclav Hapla 
5966a7c706bSVaclav Hapla    Input Parameters:
5976a7c706bSVaclav Hapla +  n  - number of values
5986a7c706bSVaclav Hapla -  X  - array of integers
5996a7c706bSVaclav Hapla 
6006a7c706bSVaclav Hapla    Output Parameters:
6016a7c706bSVaclav Hapla .  sorted - flag whether the array is sorted
6026a7c706bSVaclav Hapla 
6036a7c706bSVaclav Hapla    Level: intermediate
6046a7c706bSVaclav Hapla 
605676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortMPIInt(), PetscSortedInt(), PetscSortedReal()
6066a7c706bSVaclav Hapla @*/
6076a7c706bSVaclav Hapla PetscErrorCode  PetscSortedMPIInt(PetscInt n,const PetscMPIInt X[],PetscBool *sorted)
6086a7c706bSVaclav Hapla {
6096a7c706bSVaclav Hapla   PetscFunctionBegin;
6106a7c706bSVaclav Hapla   PetscSorted(n,X,*sorted);
6116a7c706bSVaclav Hapla   PetscFunctionReturn(0);
6126a7c706bSVaclav Hapla }
6136a7c706bSVaclav Hapla 
6146a7c706bSVaclav Hapla /*@
61517d7d925SStefano Zampini    PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order.
61617d7d925SStefano Zampini 
61717d7d925SStefano Zampini    Not Collective
61817d7d925SStefano Zampini 
61917d7d925SStefano Zampini    Input Parameters:
62017d7d925SStefano Zampini +  n  - number of values
62159e16d97SJunchao Zhang -  X  - array of integers
62217d7d925SStefano Zampini 
62317d7d925SStefano Zampini    Level: intermediate
62417d7d925SStefano Zampini 
625676f2a66SJacob Faibussowitsch    Notes:
626676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
627a5b23f4aSJose E. Roman    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
628676f2a66SJacob Faibussowitsch    code to see which routine is fastest.
629676f2a66SJacob Faibussowitsch 
630676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation()
63117d7d925SStefano Zampini @*/
632b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIInt(PetscInt n,PetscMPIInt X[])
63317d7d925SStefano Zampini {
63459e16d97SJunchao Zhang   PetscErrorCode ierr;
63559e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
63617d7d925SStefano Zampini 
63717d7d925SStefano Zampini   PetscFunctionBegin;
63859e16d97SJunchao Zhang   QuickSort1(PetscSortMPIInt,X,n,pivot,t1,ierr);
63917d7d925SStefano Zampini   PetscFunctionReturn(0);
64017d7d925SStefano Zampini }
64117d7d925SStefano Zampini 
64217d7d925SStefano Zampini /*@
64317d7d925SStefano Zampini    PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries
64417d7d925SStefano Zampini 
64517d7d925SStefano Zampini    Not Collective
64617d7d925SStefano Zampini 
64717d7d925SStefano Zampini    Input Parameters:
64817d7d925SStefano Zampini +  n  - number of values
64959e16d97SJunchao Zhang -  X  - array of integers
65017d7d925SStefano Zampini 
65117d7d925SStefano Zampini    Output Parameter:
65217d7d925SStefano Zampini .  n - number of non-redundant values
65317d7d925SStefano Zampini 
65417d7d925SStefano Zampini    Level: intermediate
65517d7d925SStefano Zampini 
65617d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
65717d7d925SStefano Zampini @*/
65859e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt X[])
65917d7d925SStefano Zampini {
66017d7d925SStefano Zampini   PetscErrorCode ierr;
66117d7d925SStefano Zampini   PetscInt       i,s = 0,N = *n, b = 0;
66217d7d925SStefano Zampini 
66317d7d925SStefano Zampini   PetscFunctionBegin;
66459e16d97SJunchao Zhang   ierr = PetscSortMPIInt(N,X);CHKERRQ(ierr);
66517d7d925SStefano Zampini   for (i=0; i<N-1; i++) {
66659e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
66759e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
668a297a907SKarl Rupp     } else s++;
66917d7d925SStefano Zampini   }
67017d7d925SStefano Zampini   *n = N - s;
67117d7d925SStefano Zampini   PetscFunctionReturn(0);
67217d7d925SStefano Zampini }
673c1f0200aSDmitry Karpeev 
6744d615ea0SBarry Smith /*@
6754d615ea0SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order;
6764d615ea0SBarry Smith        changes a second array to match the sorted first array.
6774d615ea0SBarry Smith 
6784d615ea0SBarry Smith    Not Collective
6794d615ea0SBarry Smith 
6804d615ea0SBarry Smith    Input Parameters:
6814d615ea0SBarry Smith +  n  - number of values
68259e16d97SJunchao Zhang .  X  - array of integers
68359e16d97SJunchao Zhang -  Y  - second array of integers
6844d615ea0SBarry Smith 
6854d615ea0SBarry Smith    Level: intermediate
6864d615ea0SBarry Smith 
68766f49146SPierre Jolivet .seealso: PetscMPIIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
6884d615ea0SBarry Smith @*/
689b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt X[],PetscMPIInt Y[])
6904d615ea0SBarry Smith {
6914d615ea0SBarry Smith   PetscErrorCode ierr;
69259e16d97SJunchao Zhang   PetscMPIInt    pivot,t1,t2;
6934d615ea0SBarry Smith 
6944d615ea0SBarry Smith   PetscFunctionBegin;
69559e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithArray,X,Y,n,pivot,t1,t2,ierr);
6965569a785SJunchao Zhang   PetscFunctionReturn(0);
6975569a785SJunchao Zhang }
6985569a785SJunchao Zhang 
6995569a785SJunchao Zhang /*@
7005569a785SJunchao Zhang    PetscSortMPIIntWithIntArray - Sorts an array of MPI integers in place in increasing order;
7015569a785SJunchao Zhang        changes a second array of Petsc intergers to match the sorted first array.
7025569a785SJunchao Zhang 
7035569a785SJunchao Zhang    Not Collective
7045569a785SJunchao Zhang 
7055569a785SJunchao Zhang    Input Parameters:
7065569a785SJunchao Zhang +  n  - number of values
70759e16d97SJunchao Zhang .  X  - array of MPI integers
70859e16d97SJunchao Zhang -  Y  - second array of Petsc integers
7095569a785SJunchao Zhang 
7105569a785SJunchao Zhang    Level: intermediate
7115569a785SJunchao Zhang 
7125569a785SJunchao Zhang    Notes: this routine is useful when one needs to sort MPI ranks with other integer arrays.
7135569a785SJunchao Zhang 
714676f2a66SJacob Faibussowitsch .seealso: PetscSortMPIIntWithArray(), PetscIntSortSemiOrderedWithArray(), PetscTimSortWithArray()
7155569a785SJunchao Zhang @*/
716b4f531b9SJunchao Zhang PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n,PetscMPIInt X[],PetscInt Y[])
7175569a785SJunchao Zhang {
7185569a785SJunchao Zhang   PetscErrorCode ierr;
71959e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
7205569a785SJunchao Zhang   PetscInt       t2;
7215569a785SJunchao Zhang 
7225569a785SJunchao Zhang   PetscFunctionBegin;
72359e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithIntArray,X,Y,n,pivot,t1,t2,ierr);
724e5c89e4eSSatish Balay   PetscFunctionReturn(0);
725e5c89e4eSSatish Balay }
726e5c89e4eSSatish Balay 
727e5c89e4eSSatish Balay /*@
728e5c89e4eSSatish Balay    PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
729e5c89e4eSSatish Balay        changes a second SCALAR array to match the sorted first INTEGER array.
730e5c89e4eSSatish Balay 
731e5c89e4eSSatish Balay    Not Collective
732e5c89e4eSSatish Balay 
733e5c89e4eSSatish Balay    Input Parameters:
734e5c89e4eSSatish Balay +  n  - number of values
73559e16d97SJunchao Zhang .  X  - array of integers
73659e16d97SJunchao Zhang -  Y  - second array of scalars
737e5c89e4eSSatish Balay 
738e5c89e4eSSatish Balay    Level: intermediate
739e5c89e4eSSatish Balay 
74066f49146SPierre Jolivet .seealso: PetscTimSortWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
741e5c89e4eSSatish Balay @*/
742b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithScalarArray(PetscInt n,PetscInt X[],PetscScalar Y[])
743e5c89e4eSSatish Balay {
744e5c89e4eSSatish Balay   PetscErrorCode ierr;
74559e16d97SJunchao Zhang   PetscInt       pivot,t1;
74659e16d97SJunchao Zhang   PetscScalar    t2;
747e5c89e4eSSatish Balay 
748e5c89e4eSSatish Balay   PetscFunctionBegin;
74959e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithScalarArray,X,Y,n,pivot,t1,t2,ierr);
75017df18f2SToby Isaac   PetscFunctionReturn(0);
75117df18f2SToby Isaac }
75217df18f2SToby Isaac 
7536c2863d0SBarry Smith /*@C
75417df18f2SToby Isaac    PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order;
75517df18f2SToby Isaac        changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
75617df18f2SToby Isaac        provide workspace (the size of an element in the data array) to use when sorting.
75717df18f2SToby Isaac 
75817df18f2SToby Isaac    Not Collective
75917df18f2SToby Isaac 
76017df18f2SToby Isaac    Input Parameters:
76117df18f2SToby Isaac +  n  - number of values
76259e16d97SJunchao Zhang .  X  - array of integers
76359e16d97SJunchao Zhang .  Y  - second array of data
76417df18f2SToby Isaac .  size - sizeof elements in the data array in bytes
76559e16d97SJunchao Zhang -  t2   - workspace of "size" bytes used when sorting
76617df18f2SToby Isaac 
76717df18f2SToby Isaac    Level: intermediate
76817df18f2SToby Isaac 
76966f49146SPierre Jolivet .seealso: PetscTimSortWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
77017df18f2SToby Isaac @*/
771b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithDataArray(PetscInt n,PetscInt X[],void *Y,size_t size,void *t2)
77217df18f2SToby Isaac {
77317df18f2SToby Isaac   PetscErrorCode ierr;
77459e16d97SJunchao Zhang   char           *YY = (char*)Y;
77559e16d97SJunchao Zhang   PetscInt       i,j,p,t1,pivot,hi=n-1,l,r;
77617df18f2SToby Isaac 
77717df18f2SToby Isaac   PetscFunctionBegin;
77817df18f2SToby Isaac   if (n<8) {
77959e16d97SJunchao Zhang     for (i=0; i<n; i++) {
78059e16d97SJunchao Zhang       pivot = X[i];
78159e16d97SJunchao Zhang       for (j=i+1; j<n; j++) {
78259e16d97SJunchao Zhang         if (pivot > X[j]) {
78359e16d97SJunchao Zhang           SWAP2Data(X[i],X[j],YY+size*i,YY+size*j,t1,t2,size);
78459e16d97SJunchao Zhang           pivot = X[i];
78517df18f2SToby Isaac         }
78617df18f2SToby Isaac       }
78717df18f2SToby Isaac     }
78817df18f2SToby Isaac   } else {
78959e16d97SJunchao Zhang     /* Two way partition */
79059e16d97SJunchao Zhang     p     = MEDIAN(X,hi);
79159e16d97SJunchao Zhang     pivot = X[p];
79259e16d97SJunchao Zhang     l     = 0;
79359e16d97SJunchao Zhang     r     = hi;
79459e16d97SJunchao Zhang     while (1) {
79559e16d97SJunchao Zhang       while (X[l] < pivot) l++;
79659e16d97SJunchao Zhang       while (X[r] > pivot) r--;
79759e16d97SJunchao Zhang       if (l >= r) {r++; break;}
79859e16d97SJunchao Zhang       SWAP2Data(X[l],X[r],YY+size*l,YY+size*r,t1,t2,size);
79959e16d97SJunchao Zhang       l++;
80059e16d97SJunchao Zhang       r--;
80159e16d97SJunchao Zhang     }
80259e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(l,X,Y,size,t2);CHKERRQ(ierr);
80359e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(hi-r+1,X+r,YY+size*r,size,t2);CHKERRQ(ierr);
80417df18f2SToby Isaac   }
80517df18f2SToby Isaac   PetscFunctionReturn(0);
80617df18f2SToby Isaac }
80717df18f2SToby Isaac 
80821e72a00SBarry Smith /*@
80921e72a00SBarry Smith    PetscMergeIntArray -     Merges two SORTED integer arrays, removes duplicate elements.
81021e72a00SBarry Smith 
81121e72a00SBarry Smith    Not Collective
81221e72a00SBarry Smith 
81321e72a00SBarry Smith    Input Parameters:
81421e72a00SBarry Smith +  an  - number of values in the first array
81521e72a00SBarry Smith .  aI  - first sorted array of integers
81621e72a00SBarry Smith .  bn  - number of values in the second array
81721e72a00SBarry Smith -  bI  - second array of integers
81821e72a00SBarry Smith 
81921e72a00SBarry Smith    Output Parameters:
82021e72a00SBarry Smith +  n   - number of values in the merged array
8216c2863d0SBarry Smith -  L   - merged sorted array, this is allocated if an array is not provided
82221e72a00SBarry Smith 
82321e72a00SBarry Smith    Level: intermediate
82421e72a00SBarry Smith 
82566f49146SPierre Jolivet .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
82621e72a00SBarry Smith @*/
8276c2863d0SBarry Smith PetscErrorCode  PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
82821e72a00SBarry Smith {
82921e72a00SBarry Smith   PetscErrorCode ierr;
83021e72a00SBarry Smith   PetscInt       *L_ = *L, ak, bk, k;
83121e72a00SBarry Smith 
832362febeeSStefano Zampini   PetscFunctionBegin;
83321e72a00SBarry Smith   if (!L_) {
83421e72a00SBarry Smith     ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr);
83521e72a00SBarry Smith     L_   = *L;
83621e72a00SBarry Smith   }
83721e72a00SBarry Smith   k = ak = bk = 0;
83821e72a00SBarry Smith   while (ak < an && bk < bn) {
83921e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
84021e72a00SBarry Smith       L_[k] = aI[ak];
84121e72a00SBarry Smith       ++ak;
84221e72a00SBarry Smith       ++bk;
84321e72a00SBarry Smith       ++k;
84421e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
84521e72a00SBarry Smith       L_[k] = aI[ak];
84621e72a00SBarry Smith       ++ak;
84721e72a00SBarry Smith       ++k;
84821e72a00SBarry Smith     } else {
84921e72a00SBarry Smith       L_[k] = bI[bk];
85021e72a00SBarry Smith       ++bk;
85121e72a00SBarry Smith       ++k;
85221e72a00SBarry Smith     }
85321e72a00SBarry Smith   }
85421e72a00SBarry Smith   if (ak < an) {
855580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
85621e72a00SBarry Smith     k   += (an-ak);
85721e72a00SBarry Smith   }
85821e72a00SBarry Smith   if (bk < bn) {
859580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
86021e72a00SBarry Smith     k   += (bn-bk);
86121e72a00SBarry Smith   }
86221e72a00SBarry Smith   *n = k;
86321e72a00SBarry Smith   PetscFunctionReturn(0);
86421e72a00SBarry Smith }
865b4301105SBarry Smith 
866c1f0200aSDmitry Karpeev /*@
86721e72a00SBarry Smith    PetscMergeIntArrayPair -     Merges two SORTED integer arrays that share NO common values along with an additional array of integers.
868c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
869c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
870c1f0200aSDmitry Karpeev 
871c1f0200aSDmitry Karpeev    Not Collective
872c1f0200aSDmitry Karpeev 
873c1f0200aSDmitry Karpeev    Input Parameters:
874c1f0200aSDmitry Karpeev +  an  - number of values in the first array
875c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
876c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
877c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
878c1f0200aSDmitry Karpeev .  bI  - second array of integers
879c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
880c1f0200aSDmitry Karpeev 
881c1f0200aSDmitry Karpeev    Output Parameters:
882c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
88314c49006SJed Brown .  L   - merged sorted array
884c1f0200aSDmitry Karpeev -  J   - merged additional array
885c1f0200aSDmitry Karpeev 
88695452b02SPatrick Sanan    Notes:
88795452b02SPatrick 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
888c1f0200aSDmitry Karpeev    Level: intermediate
889c1f0200aSDmitry Karpeev 
89066f49146SPierre Jolivet .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
891c1f0200aSDmitry Karpeev @*/
8926c2863d0SBarry 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)
893c1f0200aSDmitry Karpeev {
894c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
89570d8d27cSBarry Smith   PetscInt       n_, *L_, *J_, ak, bk, k;
896c1f0200aSDmitry Karpeev 
89714c49006SJed Brown   PetscFunctionBegin;
89870d8d27cSBarry Smith   PetscValidIntPointer(L,8);
89970d8d27cSBarry Smith   PetscValidIntPointer(J,9);
900c1f0200aSDmitry Karpeev   n_ = an + bn;
901c1f0200aSDmitry Karpeev   *n = n_;
90270d8d27cSBarry Smith   if (!*L) {
903785e854fSJed Brown     ierr = PetscMalloc1(n_, L);CHKERRQ(ierr);
90470d8d27cSBarry Smith   }
905d7aa01a8SHong Zhang   L_ = *L;
90670d8d27cSBarry Smith   if (!*J) {
90770d8d27cSBarry Smith     ierr = PetscMalloc1(n_, J);CHKERRQ(ierr);
908c1f0200aSDmitry Karpeev   }
909c1f0200aSDmitry Karpeev   J_   = *J;
910c1f0200aSDmitry Karpeev   k = ak = bk = 0;
911c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
912c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
913d7aa01a8SHong Zhang       L_[k] = aI[ak];
914c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
915c1f0200aSDmitry Karpeev       ++ak;
916c1f0200aSDmitry Karpeev       ++k;
917a297a907SKarl Rupp     } else {
918d7aa01a8SHong Zhang       L_[k] = bI[bk];
919c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
920c1f0200aSDmitry Karpeev       ++bk;
921c1f0200aSDmitry Karpeev       ++k;
922c1f0200aSDmitry Karpeev     }
923c1f0200aSDmitry Karpeev   }
924c1f0200aSDmitry Karpeev   if (ak < an) {
925580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
926580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,aJ+ak,an-ak);CHKERRQ(ierr);
927c1f0200aSDmitry Karpeev     k   += (an-ak);
928c1f0200aSDmitry Karpeev   }
929c1f0200aSDmitry Karpeev   if (bk < bn) {
930580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
931580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,bJ+bk,bn-bk);CHKERRQ(ierr);
932c1f0200aSDmitry Karpeev   }
933c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
934c1f0200aSDmitry Karpeev }
935c1f0200aSDmitry Karpeev 
936e498c390SJed Brown /*@
937e498c390SJed Brown    PetscMergeMPIIntArray -     Merges two SORTED integer arrays.
938e498c390SJed Brown 
939e498c390SJed Brown    Not Collective
940e498c390SJed Brown 
941e498c390SJed Brown    Input Parameters:
942e498c390SJed Brown +  an  - number of values in the first array
943e498c390SJed Brown .  aI  - first sorted array of integers
944e498c390SJed Brown .  bn  - number of values in the second array
945e498c390SJed Brown -  bI  - second array of integers
946e498c390SJed Brown 
947e498c390SJed Brown    Output Parameters:
948e498c390SJed Brown +  n   - number of values in the merged array (<= an + bn)
949e498c390SJed Brown -  L   - merged sorted array, allocated if address of NULL pointer is passed
950e498c390SJed Brown 
951e498c390SJed Brown    Level: intermediate
952e498c390SJed Brown 
95366f49146SPierre Jolivet .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
954e498c390SJed Brown @*/
955e498c390SJed Brown PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L)
956e498c390SJed Brown {
957e498c390SJed Brown   PetscErrorCode ierr;
958e498c390SJed Brown   PetscInt       ai,bi,k;
959e498c390SJed Brown 
960e498c390SJed Brown   PetscFunctionBegin;
961e498c390SJed Brown   if (!*L) {ierr = PetscMalloc1((an+bn),L);CHKERRQ(ierr);}
962e498c390SJed Brown   for (ai=0,bi=0,k=0; ai<an || bi<bn;) {
963e498c390SJed Brown     PetscInt t = -1;
964e498c390SJed Brown     for (; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
965e498c390SJed Brown     for (; bi<bn && bI[bi] == t; bi++);
966e498c390SJed Brown     for (; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
967e498c390SJed Brown     for (; ai<an && aI[ai] == t; ai++);
968e498c390SJed Brown   }
969e498c390SJed Brown   *n = k;
970e498c390SJed Brown   PetscFunctionReturn(0);
971e498c390SJed Brown }
972e5c89e4eSSatish Balay 
9736c2863d0SBarry Smith /*@C
97442eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
97542eaa7f3SBarry Smith 
97642eaa7f3SBarry Smith    Not Collective
97742eaa7f3SBarry Smith 
97842eaa7f3SBarry Smith    Input Parameters:
97942eaa7f3SBarry Smith +  n  - number of values
98042eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
98142eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
98242eaa7f3SBarry Smith 
98342eaa7f3SBarry Smith    Output Parameters:
98442eaa7f3SBarry Smith +  Nlevels - the number of levels
98542eaa7f3SBarry Smith .  Level - for each node tells its level
98642eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
98742eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
98842eaa7f3SBarry Smith -  Column - for each id tells its column index
98942eaa7f3SBarry Smith 
9906c2863d0SBarry Smith    Level: developer
99142eaa7f3SBarry Smith 
99295452b02SPatrick Sanan    Notes:
99395452b02SPatrick Sanan     This code is not currently used
99442eaa7f3SBarry Smith 
99542eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation()
99642eaa7f3SBarry Smith @*/
9977087cfbeSBarry Smith PetscErrorCode  PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column)
99842eaa7f3SBarry Smith {
99942eaa7f3SBarry Smith   PetscInt       i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column;
100042eaa7f3SBarry Smith   PetscErrorCode ierr;
1001ace3abfcSBarry Smith   PetscBool      done = PETSC_FALSE;
100242eaa7f3SBarry Smith 
100342eaa7f3SBarry Smith   PetscFunctionBegin;
10042c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mask[0],PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set");
100542eaa7f3SBarry Smith   for (i=0; i<n; i++) {
100642eaa7f3SBarry Smith     if (mask[i]) continue;
10072c71b3e2SJacob Faibussowitsch     PetscCheckFalse(parentid[i]  == i,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent");
10082c71b3e2SJacob Faibussowitsch     PetscCheckFalse(parentid[i] && mask[parentid[i]],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked");
100942eaa7f3SBarry Smith   }
101042eaa7f3SBarry Smith 
101142eaa7f3SBarry Smith   for (i=0; i<n; i++) {
101242eaa7f3SBarry Smith     if (!mask[i]) nmask++;
101342eaa7f3SBarry Smith   }
101442eaa7f3SBarry Smith 
101542eaa7f3SBarry Smith   /* determine the level in the tree of each node */
10161795a4d1SJed Brown   ierr = PetscCalloc1(n,&level);CHKERRQ(ierr);
1017a297a907SKarl Rupp 
101842eaa7f3SBarry Smith   level[0] = 1;
101942eaa7f3SBarry Smith   while (!done) {
102042eaa7f3SBarry Smith     done = PETSC_TRUE;
102142eaa7f3SBarry Smith     for (i=0; i<n; i++) {
102242eaa7f3SBarry Smith       if (mask[i]) continue;
102342eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
102442eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
102542eaa7f3SBarry Smith     }
102642eaa7f3SBarry Smith   }
102742eaa7f3SBarry Smith   for (i=0; i<n; i++) {
102842eaa7f3SBarry Smith     level[i]--;
102942eaa7f3SBarry Smith     nlevels = PetscMax(nlevels,level[i]);
103042eaa7f3SBarry Smith   }
103142eaa7f3SBarry Smith 
103242eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
10331795a4d1SJed Brown   ierr = PetscCalloc1(nlevels,&levelcnt);CHKERRQ(ierr);
103442eaa7f3SBarry Smith   for (i=0; i<n; i++) {
103542eaa7f3SBarry Smith     if (mask[i]) continue;
103642eaa7f3SBarry Smith     levelcnt[level[i]-1]++;
103742eaa7f3SBarry Smith   }
1038a297a907SKarl Rupp   for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]);
103942eaa7f3SBarry Smith 
104042eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
1041dcca6d9dSJed Brown   ierr = PetscMalloc2(levelmax,&workid,levelmax,&workparentid);CHKERRQ(ierr);
1042785e854fSJed Brown   ierr = PetscMalloc1(nmask,&idbylevel);CHKERRQ(ierr);
104342eaa7f3SBarry Smith   for (j=1; j<=nlevels;j++) {
104442eaa7f3SBarry Smith     cnt = 0;
104542eaa7f3SBarry Smith     for (i=0; i<n; i++) {
104642eaa7f3SBarry Smith       if (mask[i]) continue;
104742eaa7f3SBarry Smith       if (level[i] != j) continue;
104842eaa7f3SBarry Smith       workid[cnt]         = i;
104942eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
105042eaa7f3SBarry Smith     }
105142eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
105242eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
105342eaa7f3SBarry Smith     ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr);
105442eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
105542eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
1056580bdb30SBarry Smith     ierr  = PetscArraycpy(idbylevel+tcnt,workid,cnt);CHKERRQ(ierr);
105742eaa7f3SBarry Smith     tcnt += cnt;
105842eaa7f3SBarry Smith   }
10592c71b3e2SJacob Faibussowitsch   PetscCheckFalse(tcnt != nmask,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes");
106042eaa7f3SBarry Smith   ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr);
106142eaa7f3SBarry Smith 
106242eaa7f3SBarry Smith   /* for each node list its column */
1063785e854fSJed Brown   ierr = PetscMalloc1(n,&column);CHKERRQ(ierr);
106442eaa7f3SBarry Smith   cnt = 0;
106542eaa7f3SBarry Smith   for (j=0; j<nlevels; j++) {
106642eaa7f3SBarry Smith     for (i=0; i<levelcnt[j]; i++) {
106742eaa7f3SBarry Smith       column[idbylevel[cnt++]] = i;
106842eaa7f3SBarry Smith     }
106942eaa7f3SBarry Smith   }
107042eaa7f3SBarry Smith 
107142eaa7f3SBarry Smith   *Nlevels   = nlevels;
107242eaa7f3SBarry Smith   *Level     = level;
107342eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
107442eaa7f3SBarry Smith   *Idbylevel = idbylevel;
107542eaa7f3SBarry Smith   *Column    = column;
107642eaa7f3SBarry Smith   PetscFunctionReturn(0);
107742eaa7f3SBarry Smith }
1078ce605777SToby Isaac 
1079ce605777SToby Isaac /*@
1080ce605777SToby Isaac   PetscParallelSortedInt - Check whether an integer array, distributed over a communicator, is globally sorted.
1081ce605777SToby Isaac 
1082ce605777SToby Isaac   Collective
1083ce605777SToby Isaac 
1084ce605777SToby Isaac   Input Parameters:
1085ce605777SToby Isaac + comm - the MPI communicator
1086ce605777SToby Isaac . n - the local number of integers
1087ce605777SToby Isaac - keys - the local array of integers
1088ce605777SToby Isaac 
1089ce605777SToby Isaac   Output Parameters:
1090ce605777SToby Isaac . is_sorted - whether the array is globally sorted
1091ce605777SToby Isaac 
1092ce605777SToby Isaac   Level: developer
1093ce605777SToby Isaac 
1094ce605777SToby Isaac .seealso: PetscParallelSortInt()
1095ce605777SToby Isaac @*/
1096ce605777SToby Isaac PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
1097ce605777SToby Isaac {
1098ce605777SToby Isaac   PetscBool      sorted;
1099ce605777SToby Isaac   PetscInt       i, min, max, prevmax;
11002a1da528SToby Isaac   PetscMPIInt    rank;
1101ce605777SToby Isaac   PetscErrorCode ierr;
1102ce605777SToby Isaac 
1103ce605777SToby Isaac   PetscFunctionBegin;
1104ce605777SToby Isaac   sorted = PETSC_TRUE;
1105ce605777SToby Isaac   min    = PETSC_MAX_INT;
1106ce605777SToby Isaac   max    = PETSC_MIN_INT;
1107ce605777SToby Isaac   if (n) {
1108ce605777SToby Isaac     min = keys[0];
1109ce605777SToby Isaac     max = keys[0];
1110ce605777SToby Isaac   }
1111ce605777SToby Isaac   for (i = 1; i < n; i++) {
1112ce605777SToby Isaac     if (keys[i] < keys[i - 1]) break;
1113ce605777SToby Isaac     min = PetscMin(min,keys[i]);
1114ce605777SToby Isaac     max = PetscMax(max,keys[i]);
1115ce605777SToby Isaac   }
1116ce605777SToby Isaac   if (i < n) sorted = PETSC_FALSE;
1117ce605777SToby Isaac   prevmax = PETSC_MIN_INT;
1118ffc4695bSBarry Smith   ierr = MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(ierr);
1119ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
1120dd400576SPatrick Sanan   if (rank == 0) prevmax = PETSC_MIN_INT;
1121ce605777SToby Isaac   if (prevmax > min) sorted = PETSC_FALSE;
1122ffc4695bSBarry Smith   ierr = MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRMPI(ierr);
1123ce605777SToby Isaac   PetscFunctionReturn(0);
1124ce605777SToby Isaac }
1125