xref: /petsc/src/sys/utils/sorti.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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     t1=a; a=b; b=t1;                                                             \
29*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMemcpy(t2,c,siz));                                              \
30*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMemcpy(c,d,siz));                                               \
31*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMemcpy(d,t2,siz));                                              \
3259e16d97SJunchao Zhang   } while (0)
33e5c89e4eSSatish Balay 
34e5c89e4eSSatish Balay /*
3559e16d97SJunchao Zhang    Partition X[lo,hi] into two parts: X[lo,l) <= pivot; X[r,hi] > pivot
36e5c89e4eSSatish Balay 
3759e16d97SJunchao Zhang    Input Parameters:
3859e16d97SJunchao Zhang     + X         - array to partition
3959e16d97SJunchao Zhang     . pivot     - a pivot of X[]
4059e16d97SJunchao Zhang     . t1        - temp variable for X
4159e16d97SJunchao Zhang     - lo,hi     - lower and upper bound of the array
4259e16d97SJunchao Zhang 
4359e16d97SJunchao Zhang    Output Parameters:
4459e16d97SJunchao Zhang     + l,r       - of type PetscInt
4559e16d97SJunchao Zhang 
4659e16d97SJunchao Zhang    Notes:
4759e16d97SJunchao Zhang     The TwoWayPartition2/3 variants also partition other arrays along with X.
4859e16d97SJunchao Zhang     These arrays can have different types, so they provide their own temp t2,t3
4959e16d97SJunchao Zhang  */
5059e16d97SJunchao Zhang #define TwoWayPartition1(X,pivot,t1,lo,hi,l,r)                                   \
5159e16d97SJunchao Zhang   do {                                                                           \
5259e16d97SJunchao Zhang     l = lo;                                                                      \
5359e16d97SJunchao Zhang     r = hi;                                                                      \
5459e16d97SJunchao Zhang     while (1) {                                                                  \
5559e16d97SJunchao Zhang       while (X[l] < pivot) l++;                                                  \
5659e16d97SJunchao Zhang       while (X[r] > pivot) r--;                                                  \
5759e16d97SJunchao Zhang       if (l >= r) {r++; break;}                                                  \
5859e16d97SJunchao Zhang       SWAP1(X[l],X[r],t1);                                                       \
5959e16d97SJunchao Zhang       l++;                                                                       \
6059e16d97SJunchao Zhang       r--;                                                                       \
6159e16d97SJunchao Zhang     }                                                                            \
6259e16d97SJunchao Zhang   } while (0)
6359e16d97SJunchao Zhang 
64ce605777SToby Isaac /*
65ce605777SToby Isaac    Partition X[lo,hi] into two parts: X[lo,l) >= pivot; X[r,hi] < pivot
66ce605777SToby Isaac 
67ce605777SToby Isaac    Input Parameters:
68ce605777SToby Isaac     + X         - array to partition
69ce605777SToby Isaac     . pivot     - a pivot of X[]
70ce605777SToby Isaac     . t1        - temp variable for X
71ce605777SToby Isaac     - lo,hi     - lower and upper bound of the array
72ce605777SToby Isaac 
73ce605777SToby Isaac    Output Parameters:
74ce605777SToby Isaac     + l,r       - of type PetscInt
75ce605777SToby Isaac 
76ce605777SToby Isaac    Notes:
77ce605777SToby Isaac     The TwoWayPartition2/3 variants also partition other arrays along with X.
78ce605777SToby Isaac     These arrays can have different types, so they provide their own temp t2,t3
79ce605777SToby Isaac  */
80ce605777SToby Isaac #define TwoWayPartitionReverse1(X,pivot,t1,lo,hi,l,r)                            \
81ce605777SToby Isaac   do {                                                                           \
82ce605777SToby Isaac     l = lo;                                                                      \
83ce605777SToby Isaac     r = hi;                                                                      \
84ce605777SToby Isaac     while (1) {                                                                  \
85ce605777SToby Isaac       while (X[l] > pivot) l++;                                                  \
86ce605777SToby Isaac       while (X[r] < pivot) r--;                                                  \
87ce605777SToby Isaac       if (l >= r) {r++; break;}                                                  \
88ce605777SToby Isaac       SWAP1(X[l],X[r],t1);                                                       \
89ce605777SToby Isaac       l++;                                                                       \
90ce605777SToby Isaac       r--;                                                                       \
91ce605777SToby Isaac     }                                                                            \
92ce605777SToby Isaac   } while (0)
93ce605777SToby Isaac 
9459e16d97SJunchao Zhang #define TwoWayPartition2(X,Y,pivot,t1,t2,lo,hi,l,r)                              \
9559e16d97SJunchao Zhang   do {                                                                           \
9659e16d97SJunchao Zhang     l = lo;                                                                      \
9759e16d97SJunchao Zhang     r = hi;                                                                      \
9859e16d97SJunchao Zhang     while (1) {                                                                  \
9959e16d97SJunchao Zhang       while (X[l] < pivot) l++;                                                  \
10059e16d97SJunchao Zhang       while (X[r] > pivot) r--;                                                  \
10159e16d97SJunchao Zhang       if (l >= r) {r++; break;}                                                  \
10259e16d97SJunchao Zhang       SWAP2(X[l],X[r],Y[l],Y[r],t1,t2);                                          \
10359e16d97SJunchao Zhang       l++;                                                                       \
10459e16d97SJunchao Zhang       r--;                                                                       \
10559e16d97SJunchao Zhang     }                                                                            \
10659e16d97SJunchao Zhang   } while (0)
10759e16d97SJunchao Zhang 
10859e16d97SJunchao Zhang #define TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,lo,hi,l,r)                         \
10959e16d97SJunchao Zhang   do {                                                                           \
11059e16d97SJunchao Zhang     l = lo;                                                                      \
11159e16d97SJunchao Zhang     r = hi;                                                                      \
11259e16d97SJunchao Zhang     while (1) {                                                                  \
11359e16d97SJunchao Zhang       while (X[l] < pivot) l++;                                                  \
11459e16d97SJunchao Zhang       while (X[r] > pivot) r--;                                                  \
11559e16d97SJunchao Zhang       if (l >= r) {r++; break;}                                                  \
11659e16d97SJunchao Zhang       SWAP3(X[l],X[r],Y[l],Y[r],Z[l],Z[r],t1,t2,t3);                             \
11759e16d97SJunchao Zhang       l++;                                                                       \
11859e16d97SJunchao Zhang       r--;                                                                       \
11959e16d97SJunchao Zhang     }                                                                            \
12059e16d97SJunchao Zhang   } while (0)
12159e16d97SJunchao Zhang 
12259e16d97SJunchao Zhang /* Templates for similar functions used below */
123*5f80ce2aSJacob Faibussowitsch #define QuickSort1(FuncName,X,n,pivot,t1)                                        \
12459e16d97SJunchao Zhang   do {                                                                           \
125981bb840SJunchao Zhang     PetscCount i,j,p,l,r,hi=n-1;                                                 \
12659e16d97SJunchao Zhang     if (n < 8) {                                                                 \
12759e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
12859e16d97SJunchao Zhang         pivot = X[i];                                                            \
12959e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
13059e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
13159e16d97SJunchao Zhang             SWAP1(X[i],X[j],t1);                                                 \
13259e16d97SJunchao Zhang             pivot = X[i];                                                        \
13359e16d97SJunchao Zhang           }                                                                      \
13459e16d97SJunchao Zhang         }                                                                        \
13559e16d97SJunchao Zhang       }                                                                          \
13659e16d97SJunchao Zhang     } else {                                                                     \
13759e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
13859e16d97SJunchao Zhang       pivot = X[p];                                                              \
13959e16d97SJunchao Zhang       TwoWayPartition1(X,pivot,t1,0,hi,l,r);                                     \
140*5f80ce2aSJacob Faibussowitsch       CHKERRQ(FuncName(l,X));                                                    \
141*5f80ce2aSJacob Faibussowitsch       CHKERRQ(FuncName(hi-r+1,X+r));                                             \
14259e16d97SJunchao Zhang     }                                                                            \
14359e16d97SJunchao Zhang   } while (0)
14459e16d97SJunchao Zhang 
145ce605777SToby Isaac /* Templates for similar functions used below */
146*5f80ce2aSJacob Faibussowitsch #define QuickSortReverse1(FuncName,X,n,pivot,t1)                                 \
147ce605777SToby Isaac   do {                                                                           \
148981bb840SJunchao Zhang     PetscCount i,j,p,l,r,hi=n-1;                                                 \
149ce605777SToby Isaac     if (n < 8) {                                                                 \
150ce605777SToby Isaac       for (i=0; i<n; i++) {                                                      \
151ce605777SToby Isaac         pivot = X[i];                                                            \
152ce605777SToby Isaac         for (j=i+1; j<n; j++) {                                                  \
153ce605777SToby Isaac           if (pivot < X[j]) {                                                    \
154ce605777SToby Isaac             SWAP1(X[i],X[j],t1);                                                 \
155ce605777SToby Isaac             pivot = X[i];                                                        \
156ce605777SToby Isaac           }                                                                      \
157ce605777SToby Isaac         }                                                                        \
158ce605777SToby Isaac       }                                                                          \
159ce605777SToby Isaac     } else {                                                                     \
160ce605777SToby Isaac       p     = MEDIAN(X,hi);                                                      \
161ce605777SToby Isaac       pivot = X[p];                                                              \
162ce605777SToby Isaac       TwoWayPartitionReverse1(X,pivot,t1,0,hi,l,r);                              \
163*5f80ce2aSJacob Faibussowitsch       CHKERRQ(FuncName(l,X));                                                    \
164*5f80ce2aSJacob Faibussowitsch       CHKERRQ(FuncName(hi-r+1,X+r));                                             \
165ce605777SToby Isaac     }                                                                            \
166ce605777SToby Isaac   } while (0)
167ce605777SToby Isaac 
168*5f80ce2aSJacob Faibussowitsch #define QuickSort2(FuncName,X,Y,n,pivot,t1,t2)                                   \
16959e16d97SJunchao Zhang   do {                                                                           \
170981bb840SJunchao Zhang     PetscCount i,j,p,l,r,hi=n-1;                                                 \
17159e16d97SJunchao Zhang     if (n < 8) {                                                                 \
17259e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
17359e16d97SJunchao Zhang         pivot = X[i];                                                            \
17459e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
17559e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
17659e16d97SJunchao Zhang             SWAP2(X[i],X[j],Y[i],Y[j],t1,t2);                                    \
17759e16d97SJunchao Zhang             pivot = X[i];                                                        \
17859e16d97SJunchao Zhang           }                                                                      \
17959e16d97SJunchao Zhang         }                                                                        \
18059e16d97SJunchao Zhang       }                                                                          \
18159e16d97SJunchao Zhang     } else {                                                                     \
18259e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
18359e16d97SJunchao Zhang       pivot = X[p];                                                              \
18459e16d97SJunchao Zhang       TwoWayPartition2(X,Y,pivot,t1,t2,0,hi,l,r);                                \
185*5f80ce2aSJacob Faibussowitsch       CHKERRQ(FuncName(l,X,Y));                                                  \
186*5f80ce2aSJacob Faibussowitsch       CHKERRQ(FuncName(hi-r+1,X+r,Y+r));                                         \
18759e16d97SJunchao Zhang     }                                                                            \
18859e16d97SJunchao Zhang   } while (0)
18959e16d97SJunchao Zhang 
190*5f80ce2aSJacob Faibussowitsch #define QuickSort3(FuncName,X,Y,Z,n,pivot,t1,t2,t3)                              \
19159e16d97SJunchao Zhang   do {                                                                           \
192981bb840SJunchao Zhang     PetscCount i,j,p,l,r,hi=n-1;                                                 \
19359e16d97SJunchao Zhang     if (n < 8) {                                                                 \
19459e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
19559e16d97SJunchao Zhang         pivot = X[i];                                                            \
19659e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
19759e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
19859e16d97SJunchao Zhang             SWAP3(X[i],X[j],Y[i],Y[j],Z[i],Z[j],t1,t2,t3);                       \
19959e16d97SJunchao Zhang             pivot = X[i];                                                        \
20059e16d97SJunchao Zhang           }                                                                      \
20159e16d97SJunchao Zhang         }                                                                        \
20259e16d97SJunchao Zhang       }                                                                          \
20359e16d97SJunchao Zhang     } else {                                                                     \
20459e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
20559e16d97SJunchao Zhang       pivot = X[p];                                                              \
20659e16d97SJunchao Zhang       TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,0,hi,l,r);                           \
207*5f80ce2aSJacob Faibussowitsch       CHKERRQ(FuncName(l,X,Y,Z));                                                \
208*5f80ce2aSJacob Faibussowitsch       CHKERRQ(FuncName(hi-r+1,X+r,Y+r,Z+r));                                     \
20959e16d97SJunchao Zhang     }                                                                            \
21059e16d97SJunchao Zhang   } while (0)
211e5c89e4eSSatish Balay 
212e5c89e4eSSatish Balay /*@
2136a7c706bSVaclav Hapla    PetscSortedInt - Determines whether the array is sorted.
2146a7c706bSVaclav Hapla 
2156a7c706bSVaclav Hapla    Not Collective
2166a7c706bSVaclav Hapla 
2176a7c706bSVaclav Hapla    Input Parameters:
2186a7c706bSVaclav Hapla +  n  - number of values
2196a7c706bSVaclav Hapla -  X  - array of integers
2206a7c706bSVaclav Hapla 
2216a7c706bSVaclav Hapla    Output Parameters:
2226a7c706bSVaclav Hapla .  sorted - flag whether the array is sorted
2236a7c706bSVaclav Hapla 
2246a7c706bSVaclav Hapla    Level: intermediate
2256a7c706bSVaclav Hapla 
2266a7c706bSVaclav Hapla .seealso: PetscSortInt(), PetscSortedMPIInt(), PetscSortedReal()
2276a7c706bSVaclav Hapla @*/
2286a7c706bSVaclav Hapla PetscErrorCode  PetscSortedInt(PetscInt n,const PetscInt X[],PetscBool *sorted)
2296a7c706bSVaclav Hapla {
2306a7c706bSVaclav Hapla   PetscFunctionBegin;
231*5f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X,2);
232*5f80ce2aSJacob Faibussowitsch   PetscValidBoolPointer(sorted,3);
2336a7c706bSVaclav Hapla   PetscSorted(n,X,*sorted);
2346a7c706bSVaclav Hapla   PetscFunctionReturn(0);
2356a7c706bSVaclav Hapla }
2366a7c706bSVaclav Hapla 
2376a7c706bSVaclav Hapla /*@
238e5c89e4eSSatish Balay    PetscSortInt - Sorts an array of integers in place in increasing order.
239e5c89e4eSSatish Balay 
240e5c89e4eSSatish Balay    Not Collective
241e5c89e4eSSatish Balay 
242e5c89e4eSSatish Balay    Input Parameters:
243e5c89e4eSSatish Balay +  n  - number of values
24459e16d97SJunchao Zhang -  X  - array of integers
245e5c89e4eSSatish Balay 
246676f2a66SJacob Faibussowitsch    Notes:
247676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscIntSortSemiOrdered(), and may perform faster especially if the array
248a5b23f4aSJose E. Roman    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
249676f2a66SJacob Faibussowitsch    code to see which routine is fastest.
250676f2a66SJacob Faibussowitsch 
251e5c89e4eSSatish Balay    Level: intermediate
252e5c89e4eSSatish Balay 
253676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation()
254e5c89e4eSSatish Balay @*/
255b4f531b9SJunchao Zhang PetscErrorCode  PetscSortInt(PetscInt n,PetscInt X[])
256e5c89e4eSSatish Balay {
25759e16d97SJunchao Zhang   PetscInt pivot,t1;
258e5c89e4eSSatish Balay 
259e5c89e4eSSatish Balay   PetscFunctionBegin;
260*5f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X,2);
261*5f80ce2aSJacob Faibussowitsch   QuickSort1(PetscSortInt,X,n,pivot,t1);
262e5c89e4eSSatish Balay   PetscFunctionReturn(0);
263e5c89e4eSSatish Balay }
264e5c89e4eSSatish Balay 
2651db36a52SBarry Smith /*@
266ce605777SToby Isaac    PetscSortReverseInt - Sorts an array of integers in place in decreasing order.
267ce605777SToby Isaac 
268ce605777SToby Isaac    Not Collective
269ce605777SToby Isaac 
270ce605777SToby Isaac    Input Parameters:
271ce605777SToby Isaac +  n  - number of values
272ce605777SToby Isaac -  X  - array of integers
273ce605777SToby Isaac 
274ce605777SToby Isaac    Level: intermediate
275ce605777SToby Isaac 
276676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithPermutation()
277ce605777SToby Isaac @*/
278ce605777SToby Isaac PetscErrorCode  PetscSortReverseInt(PetscInt n,PetscInt X[])
279ce605777SToby Isaac {
280ce605777SToby Isaac   PetscInt pivot,t1;
281ce605777SToby Isaac 
282ce605777SToby Isaac   PetscFunctionBegin;
283*5f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X,2);
284*5f80ce2aSJacob Faibussowitsch   QuickSortReverse1(PetscSortReverseInt,X,n,pivot,t1);
285ce605777SToby Isaac   PetscFunctionReturn(0);
286ce605777SToby Isaac }
287ce605777SToby Isaac 
288ce605777SToby Isaac /*@
28922ab5688SLisandro Dalcin    PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted input array
29022ab5688SLisandro Dalcin 
29122ab5688SLisandro Dalcin    Not Collective
29222ab5688SLisandro Dalcin 
29322ab5688SLisandro Dalcin    Input Parameters:
29422ab5688SLisandro Dalcin +  n  - number of values
29559e16d97SJunchao Zhang -  X  - sorted array of integers
29622ab5688SLisandro Dalcin 
29722ab5688SLisandro Dalcin    Output Parameter:
29822ab5688SLisandro Dalcin .  n - number of non-redundant values
29922ab5688SLisandro Dalcin 
30022ab5688SLisandro Dalcin    Level: intermediate
30122ab5688SLisandro Dalcin 
30222ab5688SLisandro Dalcin .seealso: PetscSortInt()
30322ab5688SLisandro Dalcin @*/
30459e16d97SJunchao Zhang PetscErrorCode  PetscSortedRemoveDupsInt(PetscInt *n,PetscInt X[])
30522ab5688SLisandro Dalcin {
30622ab5688SLisandro Dalcin   PetscInt i,s = 0,N = *n, b = 0;
30722ab5688SLisandro Dalcin 
30822ab5688SLisandro Dalcin   PetscFunctionBegin;
309*5f80ce2aSJacob Faibussowitsch   PetscValidIntPointer(n,1);
310fb61b9e4SStefano Zampini   PetscCheckSorted(*n,X);
31122ab5688SLisandro Dalcin   for (i=0; i<N-1; i++) {
31259e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
31359e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
31422ab5688SLisandro Dalcin     } else s++;
31522ab5688SLisandro Dalcin   }
31622ab5688SLisandro Dalcin   *n = N - s;
31722ab5688SLisandro Dalcin   PetscFunctionReturn(0);
31822ab5688SLisandro Dalcin }
31922ab5688SLisandro Dalcin 
32022ab5688SLisandro Dalcin /*@
321b6a18d21SVaclav Hapla    PetscSortedCheckDupsInt - Checks if a sorted integer array has duplicates
322b6a18d21SVaclav Hapla 
323b6a18d21SVaclav Hapla    Not Collective
324b6a18d21SVaclav Hapla 
325b6a18d21SVaclav Hapla    Input Parameters:
326b6a18d21SVaclav Hapla +  n  - number of values
327b6a18d21SVaclav Hapla -  X  - sorted array of integers
328b6a18d21SVaclav Hapla 
329b6a18d21SVaclav Hapla    Output Parameter:
330b6a18d21SVaclav Hapla .  dups - True if the array has dups, otherwise false
331b6a18d21SVaclav Hapla 
332b6a18d21SVaclav Hapla    Level: intermediate
333b6a18d21SVaclav Hapla 
334b6a18d21SVaclav Hapla .seealso: PetscSortInt(), PetscCheckDupsInt(), PetscSortRemoveDupsInt(), PetscSortedRemoveDupsInt()
335b6a18d21SVaclav Hapla @*/
336b6a18d21SVaclav Hapla PetscErrorCode  PetscSortedCheckDupsInt(PetscInt n,const PetscInt X[],PetscBool *flg)
337b6a18d21SVaclav Hapla {
338b6a18d21SVaclav Hapla   PetscInt i;
339b6a18d21SVaclav Hapla 
340b6a18d21SVaclav Hapla   PetscFunctionBegin;
341b6a18d21SVaclav Hapla   PetscCheckSorted(n,X);
342b6a18d21SVaclav Hapla   *flg = PETSC_FALSE;
343b6a18d21SVaclav Hapla   for (i=0; i<n-1; i++) {
344b6a18d21SVaclav Hapla     if (X[i+1] == X[i]) {
345b6a18d21SVaclav Hapla       *flg = PETSC_TRUE;
346b6a18d21SVaclav Hapla       break;
347b6a18d21SVaclav Hapla     }
348b6a18d21SVaclav Hapla   }
349b6a18d21SVaclav Hapla   PetscFunctionReturn(0);
350b6a18d21SVaclav Hapla }
351b6a18d21SVaclav Hapla 
352b6a18d21SVaclav Hapla /*@
3531db36a52SBarry Smith    PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries
3541db36a52SBarry Smith 
3551db36a52SBarry Smith    Not Collective
3561db36a52SBarry Smith 
3571db36a52SBarry Smith    Input Parameters:
3581db36a52SBarry Smith +  n  - number of values
35959e16d97SJunchao Zhang -  X  - array of integers
3601db36a52SBarry Smith 
3611db36a52SBarry Smith    Output Parameter:
3621db36a52SBarry Smith .  n - number of non-redundant values
3631db36a52SBarry Smith 
3641db36a52SBarry Smith    Level: intermediate
3651db36a52SBarry Smith 
366676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortedRemoveDupsInt()
3671db36a52SBarry Smith @*/
36859e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsInt(PetscInt *n,PetscInt X[])
3691db36a52SBarry Smith {
3701db36a52SBarry Smith   PetscFunctionBegin;
371*5f80ce2aSJacob Faibussowitsch   PetscValidIntPointer(n,1);
372*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortInt(*n,X));
373*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortedRemoveDupsInt(n,X));
3741db36a52SBarry Smith   PetscFunctionReturn(0);
3751db36a52SBarry Smith }
3761db36a52SBarry Smith 
37760e03357SMatthew G Knepley /*@
378d9f1162dSJed Brown   PetscFindInt - Finds integer in a sorted array of integers
37960e03357SMatthew G Knepley 
38060e03357SMatthew G Knepley    Not Collective
38160e03357SMatthew G Knepley 
38260e03357SMatthew G Knepley    Input Parameters:
38360e03357SMatthew G Knepley +  key - the integer to locate
384d9f1162dSJed Brown .  n   - number of values in the array
38559e16d97SJunchao Zhang -  X  - array of integers
38660e03357SMatthew G Knepley 
38760e03357SMatthew G Knepley    Output Parameter:
388d9f1162dSJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
38960e03357SMatthew G Knepley 
39060e03357SMatthew G Knepley    Level: intermediate
39160e03357SMatthew G Knepley 
392676f2a66SJacob Faibussowitsch .seealso: PetscIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
39360e03357SMatthew G Knepley @*/
39459e16d97SJunchao Zhang PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
39560e03357SMatthew G Knepley {
396d9f1162dSJed Brown   PetscInt lo = 0,hi = n;
39760e03357SMatthew G Knepley 
39860e03357SMatthew G Knepley   PetscFunctionBegin;
399*5f80ce2aSJacob Faibussowitsch   PetscValidIntPointer(loc,4);
40098781d82SMatthew G Knepley   if (!n) {*loc = -1; PetscFunctionReturn(0);}
40159e16d97SJunchao Zhang   PetscValidPointer(X,3);
4026a7c706bSVaclav Hapla   PetscCheckSorted(n,X);
403d9f1162dSJed Brown   while (hi - lo > 1) {
404d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo)/2;
40559e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
406d9f1162dSJed Brown     else               lo = mid;
40760e03357SMatthew G Knepley   }
40859e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
40960e03357SMatthew G Knepley   PetscFunctionReturn(0);
41060e03357SMatthew G Knepley }
41160e03357SMatthew G Knepley 
412d2aeb606SJed Brown /*@
413f1cab4e1SJunchao Zhang   PetscCheckDupsInt - Checks if an integer array has duplicates
414f1cab4e1SJunchao Zhang 
415f1cab4e1SJunchao Zhang    Not Collective
416f1cab4e1SJunchao Zhang 
417f1cab4e1SJunchao Zhang    Input Parameters:
418f1cab4e1SJunchao Zhang +  n  - number of values in the array
419f1cab4e1SJunchao Zhang -  X  - array of integers
420f1cab4e1SJunchao Zhang 
421f1cab4e1SJunchao Zhang    Output Parameter:
422f1cab4e1SJunchao Zhang .  dups - True if the array has dups, otherwise false
423f1cab4e1SJunchao Zhang 
424f1cab4e1SJunchao Zhang    Level: intermediate
425f1cab4e1SJunchao Zhang 
426b6a18d21SVaclav Hapla .seealso: PetscSortRemoveDupsInt(), PetscSortedCheckDupsInt()
427f1cab4e1SJunchao Zhang @*/
428f1cab4e1SJunchao Zhang PetscErrorCode PetscCheckDupsInt(PetscInt n,const PetscInt X[],PetscBool *dups)
429f1cab4e1SJunchao Zhang {
430f1cab4e1SJunchao Zhang   PetscInt   i;
431f1cab4e1SJunchao Zhang   PetscHSetI ht;
432f1cab4e1SJunchao Zhang   PetscBool  missing;
433f1cab4e1SJunchao Zhang 
434f1cab4e1SJunchao Zhang   PetscFunctionBegin;
435*5f80ce2aSJacob Faibussowitsch   if (n) PetscValidIntPointer(X,2);
436*5f80ce2aSJacob Faibussowitsch   PetscValidBoolPointer(dups,3);
437f1cab4e1SJunchao Zhang   *dups = PETSC_FALSE;
438f1cab4e1SJunchao Zhang   if (n > 1) {
439*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscHSetICreate(&ht));
440*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscHSetIResize(ht,n));
441f1cab4e1SJunchao Zhang     for (i=0; i<n; i++) {
442*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscHSetIQueryAdd(ht,X[i],&missing));
443f1cab4e1SJunchao Zhang       if (!missing) {*dups = PETSC_TRUE; break;}
444f1cab4e1SJunchao Zhang     }
445*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscHSetIDestroy(&ht));
446f1cab4e1SJunchao Zhang   }
447f1cab4e1SJunchao Zhang   PetscFunctionReturn(0);
448f1cab4e1SJunchao Zhang }
449f1cab4e1SJunchao Zhang 
450f1cab4e1SJunchao Zhang /*@
451d2aeb606SJed Brown   PetscFindMPIInt - Finds MPI integer in a sorted array of integers
452d2aeb606SJed Brown 
453d2aeb606SJed Brown    Not Collective
454d2aeb606SJed Brown 
455d2aeb606SJed Brown    Input Parameters:
456d2aeb606SJed Brown +  key - the integer to locate
457d2aeb606SJed Brown .  n   - number of values in the array
45859e16d97SJunchao Zhang -  X   - array of integers
459d2aeb606SJed Brown 
460d2aeb606SJed Brown    Output Parameter:
461d2aeb606SJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
462d2aeb606SJed Brown 
463d2aeb606SJed Brown    Level: intermediate
464d2aeb606SJed Brown 
465676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
466d2aeb606SJed Brown @*/
46759e16d97SJunchao Zhang PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
468d2aeb606SJed Brown {
469d2aeb606SJed Brown   PetscInt lo = 0,hi = n;
470d2aeb606SJed Brown 
471d2aeb606SJed Brown   PetscFunctionBegin;
472d2aeb606SJed Brown   PetscValidPointer(loc,4);
473d2aeb606SJed Brown   if (!n) {*loc = -1; PetscFunctionReturn(0);}
47459e16d97SJunchao Zhang   PetscValidPointer(X,3);
4756a7c706bSVaclav Hapla   PetscCheckSorted(n,X);
476d2aeb606SJed Brown   while (hi - lo > 1) {
477d2aeb606SJed Brown     PetscInt mid = lo + (hi - lo)/2;
47859e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
479d2aeb606SJed Brown     else               lo = mid;
480d2aeb606SJed Brown   }
48159e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
482e5c89e4eSSatish Balay   PetscFunctionReturn(0);
483e5c89e4eSSatish Balay }
484e5c89e4eSSatish Balay 
485e5c89e4eSSatish Balay /*@
486e5c89e4eSSatish Balay    PetscSortIntWithArray - Sorts an array of integers in place in increasing order;
487e5c89e4eSSatish Balay        changes a second array to match the sorted first array.
488e5c89e4eSSatish Balay 
489e5c89e4eSSatish Balay    Not Collective
490e5c89e4eSSatish Balay 
491e5c89e4eSSatish Balay    Input Parameters:
492e5c89e4eSSatish Balay +  n  - number of values
49359e16d97SJunchao Zhang .  X  - array of integers
49459e16d97SJunchao Zhang -  Y  - second array of integers
495e5c89e4eSSatish Balay 
496e5c89e4eSSatish Balay    Level: intermediate
497e5c89e4eSSatish Balay 
498981bb840SJunchao Zhang .seealso: PetscIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithCountArray()
499e5c89e4eSSatish Balay @*/
500b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithArray(PetscInt n,PetscInt X[],PetscInt Y[])
501e5c89e4eSSatish Balay {
50259e16d97SJunchao Zhang   PetscInt pivot,t1,t2;
503e5c89e4eSSatish Balay 
504e5c89e4eSSatish Balay   PetscFunctionBegin;
505*5f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithArray,X,Y,n,pivot,t1,t2);
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 {
52759e16d97SJunchao Zhang   PetscInt pivot,t1,t2,t3;
528c1f0200aSDmitry Karpeev 
529c1f0200aSDmitry Karpeev   PetscFunctionBegin;
530*5f80ce2aSJacob Faibussowitsch   QuickSort3(PetscSortIntWithArrayPair,X,Y,Z,n,pivot,t1,t2,t3);
531c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
532c1f0200aSDmitry Karpeev }
533c1f0200aSDmitry Karpeev 
53417d7d925SStefano Zampini /*@
535981bb840SJunchao Zhang    PetscSortIntWithCountArray - Sorts an array of integers in place in increasing order;
536981bb840SJunchao Zhang        changes a second array to match the sorted first array.
537981bb840SJunchao Zhang 
538981bb840SJunchao Zhang    Not Collective
539981bb840SJunchao Zhang 
540981bb840SJunchao Zhang    Input Parameters:
541981bb840SJunchao Zhang +  n  - number of values
542981bb840SJunchao Zhang .  X  - array of integers
543981bb840SJunchao Zhang -  Y  - second array of PetscCounts (signed integers)
544981bb840SJunchao Zhang 
545981bb840SJunchao Zhang    Level: intermediate
546981bb840SJunchao Zhang 
547981bb840SJunchao Zhang .seealso: PetscIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
548981bb840SJunchao Zhang @*/
549981bb840SJunchao Zhang PetscErrorCode  PetscSortIntWithCountArray(PetscCount n,PetscInt X[],PetscCount Y[])
550981bb840SJunchao Zhang {
551981bb840SJunchao Zhang   PetscInt   pivot,t1;
552981bb840SJunchao Zhang   PetscCount t2;
553981bb840SJunchao Zhang 
554981bb840SJunchao Zhang   PetscFunctionBegin;
555*5f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithCountArray,X,Y,n,pivot,t1,t2);
556981bb840SJunchao Zhang   PetscFunctionReturn(0);
557981bb840SJunchao Zhang }
558981bb840SJunchao Zhang 
559981bb840SJunchao Zhang /*@
560981bb840SJunchao Zhang    PetscSortIntWithIntCountArrayPair - Sorts an array of integers in place in increasing order;
561981bb840SJunchao Zhang        changes an integer array and a PetscCount array to match the sorted first array.
562981bb840SJunchao Zhang 
563981bb840SJunchao Zhang    Not Collective
564981bb840SJunchao Zhang 
565981bb840SJunchao Zhang    Input Parameters:
566981bb840SJunchao Zhang +  n  - number of values
567981bb840SJunchao Zhang .  X  - array of integers
568981bb840SJunchao Zhang .  Y  - second array of integers (first array of the pair)
569981bb840SJunchao Zhang -  Z  - third array of PetscCounts  (second array of the pair)
570981bb840SJunchao Zhang 
571981bb840SJunchao Zhang    Level: intermediate
572981bb840SJunchao Zhang 
573981bb840SJunchao Zhang    Notes:
574981bb840SJunchao 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.
575981bb840SJunchao Zhang 
576981bb840SJunchao Zhang .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray(), PetscIntSortSemiOrdered(), PetscSortIntWithArrayPair()
577981bb840SJunchao Zhang @*/
578981bb840SJunchao Zhang PetscErrorCode  PetscSortIntWithIntCountArrayPair(PetscCount n,PetscInt X[],PetscInt Y[],PetscCount Z[])
579981bb840SJunchao Zhang {
580981bb840SJunchao Zhang   PetscInt   pivot,t1,t2; /* pivot is take from X[], so its type is still PetscInt */
581981bb840SJunchao Zhang   PetscCount t3; /* temp for Z[] */
582981bb840SJunchao Zhang 
583981bb840SJunchao Zhang   PetscFunctionBegin;
584*5f80ce2aSJacob Faibussowitsch   QuickSort3(PetscSortIntWithIntCountArrayPair,X,Y,Z,n,pivot,t1,t2,t3);
585981bb840SJunchao Zhang   PetscFunctionReturn(0);
586981bb840SJunchao Zhang }
587981bb840SJunchao Zhang 
588981bb840SJunchao Zhang /*@
5896a7c706bSVaclav Hapla    PetscSortedMPIInt - Determines whether the array is sorted.
5906a7c706bSVaclav Hapla 
5916a7c706bSVaclav Hapla    Not Collective
5926a7c706bSVaclav Hapla 
5936a7c706bSVaclav Hapla    Input Parameters:
5946a7c706bSVaclav Hapla +  n  - number of values
5956a7c706bSVaclav Hapla -  X  - array of integers
5966a7c706bSVaclav Hapla 
5976a7c706bSVaclav Hapla    Output Parameters:
5986a7c706bSVaclav Hapla .  sorted - flag whether the array is sorted
5996a7c706bSVaclav Hapla 
6006a7c706bSVaclav Hapla    Level: intermediate
6016a7c706bSVaclav Hapla 
602676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortMPIInt(), PetscSortedInt(), PetscSortedReal()
6036a7c706bSVaclav Hapla @*/
6046a7c706bSVaclav Hapla PetscErrorCode  PetscSortedMPIInt(PetscInt n,const PetscMPIInt X[],PetscBool *sorted)
6056a7c706bSVaclav Hapla {
6066a7c706bSVaclav Hapla   PetscFunctionBegin;
6076a7c706bSVaclav Hapla   PetscSorted(n,X,*sorted);
6086a7c706bSVaclav Hapla   PetscFunctionReturn(0);
6096a7c706bSVaclav Hapla }
6106a7c706bSVaclav Hapla 
6116a7c706bSVaclav Hapla /*@
61217d7d925SStefano Zampini    PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order.
61317d7d925SStefano Zampini 
61417d7d925SStefano Zampini    Not Collective
61517d7d925SStefano Zampini 
61617d7d925SStefano Zampini    Input Parameters:
61717d7d925SStefano Zampini +  n  - number of values
61859e16d97SJunchao Zhang -  X  - array of integers
61917d7d925SStefano Zampini 
62017d7d925SStefano Zampini    Level: intermediate
62117d7d925SStefano Zampini 
622676f2a66SJacob Faibussowitsch    Notes:
623676f2a66SJacob Faibussowitsch    This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
624a5b23f4aSJose E. Roman    is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
625676f2a66SJacob Faibussowitsch    code to see which routine is fastest.
626676f2a66SJacob Faibussowitsch 
627676f2a66SJacob Faibussowitsch .seealso: PetscMPIIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation()
62817d7d925SStefano Zampini @*/
629b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIInt(PetscInt n,PetscMPIInt X[])
63017d7d925SStefano Zampini {
63159e16d97SJunchao Zhang   PetscMPIInt pivot,t1;
63217d7d925SStefano Zampini 
63317d7d925SStefano Zampini   PetscFunctionBegin;
634*5f80ce2aSJacob Faibussowitsch   QuickSort1(PetscSortMPIInt,X,n,pivot,t1);
63517d7d925SStefano Zampini   PetscFunctionReturn(0);
63617d7d925SStefano Zampini }
63717d7d925SStefano Zampini 
63817d7d925SStefano Zampini /*@
63917d7d925SStefano Zampini    PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries
64017d7d925SStefano Zampini 
64117d7d925SStefano Zampini    Not Collective
64217d7d925SStefano Zampini 
64317d7d925SStefano Zampini    Input Parameters:
64417d7d925SStefano Zampini +  n  - number of values
64559e16d97SJunchao Zhang -  X  - array of integers
64617d7d925SStefano Zampini 
64717d7d925SStefano Zampini    Output Parameter:
64817d7d925SStefano Zampini .  n - number of non-redundant values
64917d7d925SStefano Zampini 
65017d7d925SStefano Zampini    Level: intermediate
65117d7d925SStefano Zampini 
65217d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
65317d7d925SStefano Zampini @*/
65459e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt X[])
65517d7d925SStefano Zampini {
656*5f80ce2aSJacob Faibussowitsch   PetscInt s = 0,N = *n,b = 0;
65717d7d925SStefano Zampini 
65817d7d925SStefano Zampini   PetscFunctionBegin;
659*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortMPIInt(N,X));
660*5f80ce2aSJacob Faibussowitsch   for (PetscInt i=0; i<N-1; i++) {
66159e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
66259e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
663a297a907SKarl Rupp     } else s++;
66417d7d925SStefano Zampini   }
66517d7d925SStefano Zampini   *n = N - s;
66617d7d925SStefano Zampini   PetscFunctionReturn(0);
66717d7d925SStefano Zampini }
668c1f0200aSDmitry Karpeev 
6694d615ea0SBarry Smith /*@
6704d615ea0SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order;
6714d615ea0SBarry Smith        changes a second array to match the sorted first array.
6724d615ea0SBarry Smith 
6734d615ea0SBarry Smith    Not Collective
6744d615ea0SBarry Smith 
6754d615ea0SBarry Smith    Input Parameters:
6764d615ea0SBarry Smith +  n  - number of values
67759e16d97SJunchao Zhang .  X  - array of integers
67859e16d97SJunchao Zhang -  Y  - second array of integers
6794d615ea0SBarry Smith 
6804d615ea0SBarry Smith    Level: intermediate
6814d615ea0SBarry Smith 
68266f49146SPierre Jolivet .seealso: PetscMPIIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
6834d615ea0SBarry Smith @*/
684b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt X[],PetscMPIInt Y[])
6854d615ea0SBarry Smith {
68659e16d97SJunchao Zhang   PetscMPIInt pivot,t1,t2;
6874d615ea0SBarry Smith 
6884d615ea0SBarry Smith   PetscFunctionBegin;
689*5f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortMPIIntWithArray,X,Y,n,pivot,t1,t2);
6905569a785SJunchao Zhang   PetscFunctionReturn(0);
6915569a785SJunchao Zhang }
6925569a785SJunchao Zhang 
6935569a785SJunchao Zhang /*@
6945569a785SJunchao Zhang    PetscSortMPIIntWithIntArray - Sorts an array of MPI integers in place in increasing order;
6955569a785SJunchao Zhang        changes a second array of Petsc intergers to match the sorted first array.
6965569a785SJunchao Zhang 
6975569a785SJunchao Zhang    Not Collective
6985569a785SJunchao Zhang 
6995569a785SJunchao Zhang    Input Parameters:
7005569a785SJunchao Zhang +  n  - number of values
70159e16d97SJunchao Zhang .  X  - array of MPI integers
70259e16d97SJunchao Zhang -  Y  - second array of Petsc integers
7035569a785SJunchao Zhang 
7045569a785SJunchao Zhang    Level: intermediate
7055569a785SJunchao Zhang 
7065569a785SJunchao Zhang    Notes: this routine is useful when one needs to sort MPI ranks with other integer arrays.
7075569a785SJunchao Zhang 
708676f2a66SJacob Faibussowitsch .seealso: PetscSortMPIIntWithArray(), PetscIntSortSemiOrderedWithArray(), PetscTimSortWithArray()
7095569a785SJunchao Zhang @*/
710b4f531b9SJunchao Zhang PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n,PetscMPIInt X[],PetscInt Y[])
7115569a785SJunchao Zhang {
71259e16d97SJunchao Zhang   PetscMPIInt pivot,t1;
7135569a785SJunchao Zhang   PetscInt    t2;
7145569a785SJunchao Zhang 
7155569a785SJunchao Zhang   PetscFunctionBegin;
716*5f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortMPIIntWithIntArray,X,Y,n,pivot,t1,t2);
717e5c89e4eSSatish Balay   PetscFunctionReturn(0);
718e5c89e4eSSatish Balay }
719e5c89e4eSSatish Balay 
720e5c89e4eSSatish Balay /*@
721e5c89e4eSSatish Balay    PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
722e5c89e4eSSatish Balay        changes a second SCALAR array to match the sorted first INTEGER array.
723e5c89e4eSSatish Balay 
724e5c89e4eSSatish Balay    Not Collective
725e5c89e4eSSatish Balay 
726e5c89e4eSSatish Balay    Input Parameters:
727e5c89e4eSSatish Balay +  n  - number of values
72859e16d97SJunchao Zhang .  X  - array of integers
72959e16d97SJunchao Zhang -  Y  - second array of scalars
730e5c89e4eSSatish Balay 
731e5c89e4eSSatish Balay    Level: intermediate
732e5c89e4eSSatish Balay 
73366f49146SPierre Jolivet .seealso: PetscTimSortWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
734e5c89e4eSSatish Balay @*/
735b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithScalarArray(PetscInt n,PetscInt X[],PetscScalar Y[])
736e5c89e4eSSatish Balay {
73759e16d97SJunchao Zhang   PetscInt    pivot,t1;
73859e16d97SJunchao Zhang   PetscScalar t2;
739e5c89e4eSSatish Balay 
740e5c89e4eSSatish Balay   PetscFunctionBegin;
741*5f80ce2aSJacob Faibussowitsch   QuickSort2(PetscSortIntWithScalarArray,X,Y,n,pivot,t1,t2);
74217df18f2SToby Isaac   PetscFunctionReturn(0);
74317df18f2SToby Isaac }
74417df18f2SToby Isaac 
7456c2863d0SBarry Smith /*@C
74617df18f2SToby Isaac    PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order;
74717df18f2SToby Isaac        changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
74817df18f2SToby Isaac        provide workspace (the size of an element in the data array) to use when sorting.
74917df18f2SToby Isaac 
75017df18f2SToby Isaac    Not Collective
75117df18f2SToby Isaac 
75217df18f2SToby Isaac    Input Parameters:
75317df18f2SToby Isaac +  n  - number of values
75459e16d97SJunchao Zhang .  X  - array of integers
75559e16d97SJunchao Zhang .  Y  - second array of data
75617df18f2SToby Isaac .  size - sizeof elements in the data array in bytes
75759e16d97SJunchao Zhang -  t2   - workspace of "size" bytes used when sorting
75817df18f2SToby Isaac 
75917df18f2SToby Isaac    Level: intermediate
76017df18f2SToby Isaac 
76166f49146SPierre Jolivet .seealso: PetscTimSortWithArray(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
76217df18f2SToby Isaac @*/
763b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithDataArray(PetscInt n,PetscInt X[],void *Y,size_t size,void *t2)
76417df18f2SToby Isaac {
76559e16d97SJunchao Zhang   char     *YY       = (char*)Y;
766*5f80ce2aSJacob Faibussowitsch   PetscInt  t1,pivot,hi = n-1;
76717df18f2SToby Isaac 
76817df18f2SToby Isaac   PetscFunctionBegin;
76917df18f2SToby Isaac   if (n<8) {
770*5f80ce2aSJacob Faibussowitsch     for (PetscInt i=0; i<n; i++) {
77159e16d97SJunchao Zhang       pivot = X[i];
772*5f80ce2aSJacob Faibussowitsch       for (PetscInt j=i+1; j<n; j++) {
77359e16d97SJunchao Zhang         if (pivot > X[j]) {
77459e16d97SJunchao Zhang           SWAP2Data(X[i],X[j],YY+size*i,YY+size*j,t1,t2,size);
77559e16d97SJunchao Zhang           pivot = X[i];
77617df18f2SToby Isaac         }
77717df18f2SToby Isaac       }
77817df18f2SToby Isaac     }
77917df18f2SToby Isaac   } else {
78059e16d97SJunchao Zhang     /* Two way partition */
781*5f80ce2aSJacob Faibussowitsch     PetscInt l = 0,r = hi;
782*5f80ce2aSJacob Faibussowitsch 
783*5f80ce2aSJacob Faibussowitsch     pivot = X[MEDIAN(X,hi)];
78459e16d97SJunchao Zhang     while (1) {
78559e16d97SJunchao Zhang       while (X[l] < pivot) l++;
78659e16d97SJunchao Zhang       while (X[r] > pivot) r--;
78759e16d97SJunchao Zhang       if (l >= r) {r++; break;}
78859e16d97SJunchao Zhang       SWAP2Data(X[l],X[r],YY+size*l,YY+size*r,t1,t2,size);
78959e16d97SJunchao Zhang       l++;
79059e16d97SJunchao Zhang       r--;
79159e16d97SJunchao Zhang     }
792*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSortIntWithDataArray(l,X,Y,size,t2));
793*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSortIntWithDataArray(hi-r+1,X+r,YY+size*r,size,t2));
79417df18f2SToby Isaac   }
79517df18f2SToby Isaac   PetscFunctionReturn(0);
79617df18f2SToby Isaac }
79717df18f2SToby Isaac 
79821e72a00SBarry Smith /*@
79921e72a00SBarry Smith    PetscMergeIntArray -     Merges two SORTED integer arrays, removes duplicate elements.
80021e72a00SBarry Smith 
80121e72a00SBarry Smith    Not Collective
80221e72a00SBarry Smith 
80321e72a00SBarry Smith    Input Parameters:
80421e72a00SBarry Smith +  an  - number of values in the first array
80521e72a00SBarry Smith .  aI  - first sorted array of integers
80621e72a00SBarry Smith .  bn  - number of values in the second array
80721e72a00SBarry Smith -  bI  - second array of integers
80821e72a00SBarry Smith 
80921e72a00SBarry Smith    Output Parameters:
81021e72a00SBarry Smith +  n   - number of values in the merged array
8116c2863d0SBarry Smith -  L   - merged sorted array, this is allocated if an array is not provided
81221e72a00SBarry Smith 
81321e72a00SBarry Smith    Level: intermediate
81421e72a00SBarry Smith 
81566f49146SPierre Jolivet .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
81621e72a00SBarry Smith @*/
8176c2863d0SBarry Smith PetscErrorCode  PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
81821e72a00SBarry Smith {
81921e72a00SBarry Smith   PetscInt       *L_ = *L, ak, bk, k;
82021e72a00SBarry Smith 
821362febeeSStefano Zampini   PetscFunctionBegin;
82221e72a00SBarry Smith   if (!L_) {
823*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(an+bn, L));
82421e72a00SBarry Smith     L_   = *L;
82521e72a00SBarry Smith   }
82621e72a00SBarry Smith   k = ak = bk = 0;
82721e72a00SBarry Smith   while (ak < an && bk < bn) {
82821e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
82921e72a00SBarry Smith       L_[k] = aI[ak];
83021e72a00SBarry Smith       ++ak;
83121e72a00SBarry Smith       ++bk;
83221e72a00SBarry Smith       ++k;
83321e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
83421e72a00SBarry Smith       L_[k] = aI[ak];
83521e72a00SBarry Smith       ++ak;
83621e72a00SBarry Smith       ++k;
83721e72a00SBarry Smith     } else {
83821e72a00SBarry Smith       L_[k] = bI[bk];
83921e72a00SBarry Smith       ++bk;
84021e72a00SBarry Smith       ++k;
84121e72a00SBarry Smith     }
84221e72a00SBarry Smith   }
84321e72a00SBarry Smith   if (ak < an) {
844*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(L_+k,aI+ak,an-ak));
84521e72a00SBarry Smith     k   += (an-ak);
84621e72a00SBarry Smith   }
84721e72a00SBarry Smith   if (bk < bn) {
848*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(L_+k,bI+bk,bn-bk));
84921e72a00SBarry Smith     k   += (bn-bk);
85021e72a00SBarry Smith   }
85121e72a00SBarry Smith   *n = k;
85221e72a00SBarry Smith   PetscFunctionReturn(0);
85321e72a00SBarry Smith }
854b4301105SBarry Smith 
855c1f0200aSDmitry Karpeev /*@
85621e72a00SBarry Smith    PetscMergeIntArrayPair -     Merges two SORTED integer arrays that share NO common values along with an additional array of integers.
857c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
858c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
859c1f0200aSDmitry Karpeev 
860c1f0200aSDmitry Karpeev    Not Collective
861c1f0200aSDmitry Karpeev 
862c1f0200aSDmitry Karpeev    Input Parameters:
863c1f0200aSDmitry Karpeev +  an  - number of values in the first array
864c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
865c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
866c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
867c1f0200aSDmitry Karpeev .  bI  - second array of integers
868c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
869c1f0200aSDmitry Karpeev 
870c1f0200aSDmitry Karpeev    Output Parameters:
871c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
87214c49006SJed Brown .  L   - merged sorted array
873c1f0200aSDmitry Karpeev -  J   - merged additional array
874c1f0200aSDmitry Karpeev 
87595452b02SPatrick Sanan    Notes:
87695452b02SPatrick 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
877c1f0200aSDmitry Karpeev    Level: intermediate
878c1f0200aSDmitry Karpeev 
87966f49146SPierre Jolivet .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
880c1f0200aSDmitry Karpeev @*/
8816c2863d0SBarry 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)
882c1f0200aSDmitry Karpeev {
88370d8d27cSBarry Smith   PetscInt       n_, *L_, *J_, ak, bk, k;
884c1f0200aSDmitry Karpeev 
88514c49006SJed Brown   PetscFunctionBegin;
88670d8d27cSBarry Smith   PetscValidIntPointer(L,8);
88770d8d27cSBarry Smith   PetscValidIntPointer(J,9);
888c1f0200aSDmitry Karpeev   n_ = an + bn;
889c1f0200aSDmitry Karpeev   *n = n_;
89070d8d27cSBarry Smith   if (!*L) {
891*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(n_, L));
89270d8d27cSBarry Smith   }
893d7aa01a8SHong Zhang   L_ = *L;
89470d8d27cSBarry Smith   if (!*J) {
895*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(n_, J));
896c1f0200aSDmitry Karpeev   }
897c1f0200aSDmitry Karpeev   J_   = *J;
898c1f0200aSDmitry Karpeev   k = ak = bk = 0;
899c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
900c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
901d7aa01a8SHong Zhang       L_[k] = aI[ak];
902c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
903c1f0200aSDmitry Karpeev       ++ak;
904c1f0200aSDmitry Karpeev       ++k;
905a297a907SKarl Rupp     } else {
906d7aa01a8SHong Zhang       L_[k] = bI[bk];
907c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
908c1f0200aSDmitry Karpeev       ++bk;
909c1f0200aSDmitry Karpeev       ++k;
910c1f0200aSDmitry Karpeev     }
911c1f0200aSDmitry Karpeev   }
912c1f0200aSDmitry Karpeev   if (ak < an) {
913*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(L_+k,aI+ak,an-ak));
914*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(J_+k,aJ+ak,an-ak));
915c1f0200aSDmitry Karpeev     k   += (an-ak);
916c1f0200aSDmitry Karpeev   }
917c1f0200aSDmitry Karpeev   if (bk < bn) {
918*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(L_+k,bI+bk,bn-bk));
919*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(J_+k,bJ+bk,bn-bk));
920c1f0200aSDmitry Karpeev   }
921c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
922c1f0200aSDmitry Karpeev }
923c1f0200aSDmitry Karpeev 
924e498c390SJed Brown /*@
925e498c390SJed Brown    PetscMergeMPIIntArray -     Merges two SORTED integer arrays.
926e498c390SJed Brown 
927e498c390SJed Brown    Not Collective
928e498c390SJed Brown 
929e498c390SJed Brown    Input Parameters:
930e498c390SJed Brown +  an  - number of values in the first array
931e498c390SJed Brown .  aI  - first sorted array of integers
932e498c390SJed Brown .  bn  - number of values in the second array
933e498c390SJed Brown -  bI  - second array of integers
934e498c390SJed Brown 
935e498c390SJed Brown    Output Parameters:
936e498c390SJed Brown +  n   - number of values in the merged array (<= an + bn)
937e498c390SJed Brown -  L   - merged sorted array, allocated if address of NULL pointer is passed
938e498c390SJed Brown 
939e498c390SJed Brown    Level: intermediate
940e498c390SJed Brown 
94166f49146SPierre Jolivet .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortIntWithArray()
942e498c390SJed Brown @*/
943e498c390SJed Brown PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L)
944e498c390SJed Brown {
945e498c390SJed Brown   PetscInt       ai,bi,k;
946e498c390SJed Brown 
947e498c390SJed Brown   PetscFunctionBegin;
948*5f80ce2aSJacob Faibussowitsch   if (!*L) CHKERRQ(PetscMalloc1((an+bn),L));
949e498c390SJed Brown   for (ai=0,bi=0,k=0; ai<an || bi<bn;) {
950e498c390SJed Brown     PetscInt t = -1;
951e498c390SJed Brown     for (; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
952e498c390SJed Brown     for (; bi<bn && bI[bi] == t; bi++);
953e498c390SJed Brown     for (; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
954e498c390SJed Brown     for (; ai<an && aI[ai] == t; ai++);
955e498c390SJed Brown   }
956e498c390SJed Brown   *n = k;
957e498c390SJed Brown   PetscFunctionReturn(0);
958e498c390SJed Brown }
959e5c89e4eSSatish Balay 
9606c2863d0SBarry Smith /*@C
96142eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
96242eaa7f3SBarry Smith 
96342eaa7f3SBarry Smith    Not Collective
96442eaa7f3SBarry Smith 
96542eaa7f3SBarry Smith    Input Parameters:
96642eaa7f3SBarry Smith +  n  - number of values
96742eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
96842eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
96942eaa7f3SBarry Smith 
97042eaa7f3SBarry Smith    Output Parameters:
97142eaa7f3SBarry Smith +  Nlevels - the number of levels
97242eaa7f3SBarry Smith .  Level - for each node tells its level
97342eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
97442eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
97542eaa7f3SBarry Smith -  Column - for each id tells its column index
97642eaa7f3SBarry Smith 
9776c2863d0SBarry Smith    Level: developer
97842eaa7f3SBarry Smith 
97995452b02SPatrick Sanan    Notes:
98095452b02SPatrick Sanan     This code is not currently used
98142eaa7f3SBarry Smith 
98242eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation()
98342eaa7f3SBarry Smith @*/
9847087cfbeSBarry Smith PetscErrorCode  PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column)
98542eaa7f3SBarry Smith {
98642eaa7f3SBarry Smith   PetscInt       i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column;
987ace3abfcSBarry Smith   PetscBool      done = PETSC_FALSE;
98842eaa7f3SBarry Smith 
98942eaa7f3SBarry Smith   PetscFunctionBegin;
9902c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mask[0],PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set");
99142eaa7f3SBarry Smith   for (i=0; i<n; i++) {
99242eaa7f3SBarry Smith     if (mask[i]) continue;
9932c71b3e2SJacob Faibussowitsch     PetscCheckFalse(parentid[i]  == i,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent");
9942c71b3e2SJacob Faibussowitsch     PetscCheckFalse(parentid[i] && mask[parentid[i]],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked");
99542eaa7f3SBarry Smith   }
99642eaa7f3SBarry Smith 
99742eaa7f3SBarry Smith   for (i=0; i<n; i++) {
99842eaa7f3SBarry Smith     if (!mask[i]) nmask++;
99942eaa7f3SBarry Smith   }
100042eaa7f3SBarry Smith 
100142eaa7f3SBarry Smith   /* determine the level in the tree of each node */
1002*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(n,&level));
1003a297a907SKarl Rupp 
100442eaa7f3SBarry Smith   level[0] = 1;
100542eaa7f3SBarry Smith   while (!done) {
100642eaa7f3SBarry Smith     done = PETSC_TRUE;
100742eaa7f3SBarry Smith     for (i=0; i<n; i++) {
100842eaa7f3SBarry Smith       if (mask[i]) continue;
100942eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
101042eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
101142eaa7f3SBarry Smith     }
101242eaa7f3SBarry Smith   }
101342eaa7f3SBarry Smith   for (i=0; i<n; i++) {
101442eaa7f3SBarry Smith     level[i]--;
101542eaa7f3SBarry Smith     nlevels = PetscMax(nlevels,level[i]);
101642eaa7f3SBarry Smith   }
101742eaa7f3SBarry Smith 
101842eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
1019*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(nlevels,&levelcnt));
102042eaa7f3SBarry Smith   for (i=0; i<n; i++) {
102142eaa7f3SBarry Smith     if (mask[i]) continue;
102242eaa7f3SBarry Smith     levelcnt[level[i]-1]++;
102342eaa7f3SBarry Smith   }
1024a297a907SKarl Rupp   for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]);
102542eaa7f3SBarry Smith 
102642eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
1027*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(levelmax,&workid,levelmax,&workparentid));
1028*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nmask,&idbylevel));
102942eaa7f3SBarry Smith   for (j=1; j<=nlevels;j++) {
103042eaa7f3SBarry Smith     cnt = 0;
103142eaa7f3SBarry Smith     for (i=0; i<n; i++) {
103242eaa7f3SBarry Smith       if (mask[i]) continue;
103342eaa7f3SBarry Smith       if (level[i] != j) continue;
103442eaa7f3SBarry Smith       workid[cnt]         = i;
103542eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
103642eaa7f3SBarry Smith     }
103742eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
103842eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
1039*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSortIntWithArray(cnt,workparentid,workid));
104042eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
104142eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
1042*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(idbylevel+tcnt,workid,cnt));
104342eaa7f3SBarry Smith     tcnt += cnt;
104442eaa7f3SBarry Smith   }
10452c71b3e2SJacob Faibussowitsch   PetscCheckFalse(tcnt != nmask,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes");
1046*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(workid,workparentid));
104742eaa7f3SBarry Smith 
104842eaa7f3SBarry Smith   /* for each node list its column */
1049*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&column));
105042eaa7f3SBarry Smith   cnt = 0;
105142eaa7f3SBarry Smith   for (j=0; j<nlevels; j++) {
105242eaa7f3SBarry Smith     for (i=0; i<levelcnt[j]; i++) {
105342eaa7f3SBarry Smith       column[idbylevel[cnt++]] = i;
105442eaa7f3SBarry Smith     }
105542eaa7f3SBarry Smith   }
105642eaa7f3SBarry Smith 
105742eaa7f3SBarry Smith   *Nlevels   = nlevels;
105842eaa7f3SBarry Smith   *Level     = level;
105942eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
106042eaa7f3SBarry Smith   *Idbylevel = idbylevel;
106142eaa7f3SBarry Smith   *Column    = column;
106242eaa7f3SBarry Smith   PetscFunctionReturn(0);
106342eaa7f3SBarry Smith }
1064ce605777SToby Isaac 
1065ce605777SToby Isaac /*@
1066ce605777SToby Isaac   PetscParallelSortedInt - Check whether an integer array, distributed over a communicator, is globally sorted.
1067ce605777SToby Isaac 
1068ce605777SToby Isaac   Collective
1069ce605777SToby Isaac 
1070ce605777SToby Isaac   Input Parameters:
1071ce605777SToby Isaac + comm - the MPI communicator
1072ce605777SToby Isaac . n - the local number of integers
1073ce605777SToby Isaac - keys - the local array of integers
1074ce605777SToby Isaac 
1075ce605777SToby Isaac   Output Parameters:
1076ce605777SToby Isaac . is_sorted - whether the array is globally sorted
1077ce605777SToby Isaac 
1078ce605777SToby Isaac   Level: developer
1079ce605777SToby Isaac 
1080ce605777SToby Isaac .seealso: PetscParallelSortInt()
1081ce605777SToby Isaac @*/
1082ce605777SToby Isaac PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
1083ce605777SToby Isaac {
1084ce605777SToby Isaac   PetscBool      sorted;
1085ce605777SToby Isaac   PetscInt       i, min, max, prevmax;
10862a1da528SToby Isaac   PetscMPIInt    rank;
1087ce605777SToby Isaac 
1088ce605777SToby Isaac   PetscFunctionBegin;
1089ce605777SToby Isaac   sorted = PETSC_TRUE;
1090ce605777SToby Isaac   min    = PETSC_MAX_INT;
1091ce605777SToby Isaac   max    = PETSC_MIN_INT;
1092ce605777SToby Isaac   if (n) {
1093ce605777SToby Isaac     min = keys[0];
1094ce605777SToby Isaac     max = keys[0];
1095ce605777SToby Isaac   }
1096ce605777SToby Isaac   for (i = 1; i < n; i++) {
1097ce605777SToby Isaac     if (keys[i] < keys[i - 1]) break;
1098ce605777SToby Isaac     min = PetscMin(min,keys[i]);
1099ce605777SToby Isaac     max = PetscMax(max,keys[i]);
1100ce605777SToby Isaac   }
1101ce605777SToby Isaac   if (i < n) sorted = PETSC_FALSE;
1102ce605777SToby Isaac   prevmax = PETSC_MIN_INT;
1103*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm));
1104*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
1105dd400576SPatrick Sanan   if (rank == 0) prevmax = PETSC_MIN_INT;
1106ce605777SToby Isaac   if (prevmax > min) sorted = PETSC_FALSE;
1107*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm));
1108ce605777SToby Isaac   PetscFunctionReturn(0);
1109ce605777SToby Isaac }
1110