xref: /petsc/src/sys/utils/sorti.c (revision f1cab4e1e924af9242bb4b99cbdea45bd01822e9)
17d0a6c19SBarry Smith 
2e5c89e4eSSatish Balay /*
3e5c89e4eSSatish Balay    This file contains routines for sorting integers. Values are sorted in place.
459e16d97SJunchao Zhang    One can use src/sys/examples/tests/ex52.c for benchmarking.
5e5c89e4eSSatish Balay  */
6af0996ceSBarry Smith #include <petsc/private/petscimpl.h>                /*I  "petscsys.h"  I*/
7*f1cab4e1SJunchao 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 
6559e16d97SJunchao Zhang #define TwoWayPartition2(X,Y,pivot,t1,t2,lo,hi,l,r)                              \
6659e16d97SJunchao Zhang   do {                                                                           \
6759e16d97SJunchao Zhang     l = lo;                                                                      \
6859e16d97SJunchao Zhang     r = hi;                                                                      \
6959e16d97SJunchao Zhang     while(1) {                                                                   \
7059e16d97SJunchao Zhang       while (X[l] < pivot) l++;                                                  \
7159e16d97SJunchao Zhang       while (X[r] > pivot) r--;                                                  \
7259e16d97SJunchao Zhang       if (l >= r) {r++; break;}                                                  \
7359e16d97SJunchao Zhang       SWAP2(X[l],X[r],Y[l],Y[r],t1,t2);                                          \
7459e16d97SJunchao Zhang       l++;                                                                       \
7559e16d97SJunchao Zhang       r--;                                                                       \
7659e16d97SJunchao Zhang     }                                                                            \
7759e16d97SJunchao Zhang   } while(0)
7859e16d97SJunchao Zhang 
7959e16d97SJunchao Zhang #define TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,lo,hi,l,r)                         \
8059e16d97SJunchao Zhang   do {                                                                           \
8159e16d97SJunchao Zhang     l = lo;                                                                      \
8259e16d97SJunchao Zhang     r = hi;                                                                      \
8359e16d97SJunchao Zhang     while(1) {                                                                   \
8459e16d97SJunchao Zhang       while (X[l] < pivot) l++;                                                  \
8559e16d97SJunchao Zhang       while (X[r] > pivot) r--;                                                  \
8659e16d97SJunchao Zhang       if (l >= r) {r++; break;}                                                  \
8759e16d97SJunchao Zhang       SWAP3(X[l],X[r],Y[l],Y[r],Z[l],Z[r],t1,t2,t3);                             \
8859e16d97SJunchao Zhang       l++;                                                                       \
8959e16d97SJunchao Zhang       r--;                                                                       \
9059e16d97SJunchao Zhang     }                                                                            \
9159e16d97SJunchao Zhang   } while(0)
9259e16d97SJunchao Zhang 
9359e16d97SJunchao Zhang /* Templates for similar functions used below */
9459e16d97SJunchao Zhang #define QuickSort1(FuncName,X,n,pivot,t1,ierr)                                   \
9559e16d97SJunchao Zhang   do {                                                                           \
9659e16d97SJunchao Zhang     PetscInt i,j,p,l,r,hi=n-1;                                                   \
9759e16d97SJunchao Zhang     if (n < 8) {                                                                 \
9859e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
9959e16d97SJunchao Zhang         pivot = X[i];                                                            \
10059e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
10159e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
10259e16d97SJunchao Zhang             SWAP1(X[i],X[j],t1);                                                 \
10359e16d97SJunchao Zhang             pivot = X[i];                                                        \
10459e16d97SJunchao Zhang           }                                                                      \
10559e16d97SJunchao Zhang         }                                                                        \
10659e16d97SJunchao Zhang       }                                                                          \
10759e16d97SJunchao Zhang     } else {                                                                     \
10859e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
10959e16d97SJunchao Zhang       pivot = X[p];                                                              \
11059e16d97SJunchao Zhang       TwoWayPartition1(X,pivot,t1,0,hi,l,r);                                     \
11159e16d97SJunchao Zhang       ierr  = FuncName(l,X);CHKERRQ(ierr);                                       \
11259e16d97SJunchao Zhang       ierr  = FuncName(hi-r+1,X+r);CHKERRQ(ierr);                                \
11359e16d97SJunchao Zhang     }                                                                            \
11459e16d97SJunchao Zhang   } while(0)
11559e16d97SJunchao Zhang 
11659e16d97SJunchao Zhang #define QuickSort2(FuncName,X,Y,n,pivot,t1,t2,ierr)                              \
11759e16d97SJunchao Zhang   do {                                                                           \
11859e16d97SJunchao Zhang     PetscInt i,j,p,l,r,hi=n-1;                                                   \
11959e16d97SJunchao Zhang     if (n < 8) {                                                                 \
12059e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
12159e16d97SJunchao Zhang         pivot = X[i];                                                            \
12259e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
12359e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
12459e16d97SJunchao Zhang             SWAP2(X[i],X[j],Y[i],Y[j],t1,t2);                                    \
12559e16d97SJunchao Zhang             pivot = X[i];                                                        \
12659e16d97SJunchao Zhang           }                                                                      \
12759e16d97SJunchao Zhang         }                                                                        \
12859e16d97SJunchao Zhang       }                                                                          \
12959e16d97SJunchao Zhang     } else {                                                                     \
13059e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
13159e16d97SJunchao Zhang       pivot = X[p];                                                              \
13259e16d97SJunchao Zhang       TwoWayPartition2(X,Y,pivot,t1,t2,0,hi,l,r);                                \
13359e16d97SJunchao Zhang       ierr  = FuncName(l,X,Y);CHKERRQ(ierr);                                     \
13459e16d97SJunchao Zhang       ierr  = FuncName(hi-r+1,X+r,Y+r);CHKERRQ(ierr);                            \
13559e16d97SJunchao Zhang     }                                                                            \
13659e16d97SJunchao Zhang   } while(0)
13759e16d97SJunchao Zhang 
13859e16d97SJunchao Zhang #define QuickSort3(FuncName,X,Y,Z,n,pivot,t1,t2,t3,ierr)                         \
13959e16d97SJunchao Zhang   do {                                                                           \
14059e16d97SJunchao Zhang     PetscInt i,j,p,l,r,hi=n-1;                                                   \
14159e16d97SJunchao Zhang     if (n < 8) {                                                                 \
14259e16d97SJunchao Zhang       for (i=0; i<n; i++) {                                                      \
14359e16d97SJunchao Zhang         pivot = X[i];                                                            \
14459e16d97SJunchao Zhang         for (j=i+1; j<n; j++) {                                                  \
14559e16d97SJunchao Zhang           if (pivot > X[j]) {                                                    \
14659e16d97SJunchao Zhang             SWAP3(X[i],X[j],Y[i],Y[j],Z[i],Z[j],t1,t2,t3);                       \
14759e16d97SJunchao Zhang             pivot = X[i];                                                        \
14859e16d97SJunchao Zhang           }                                                                      \
14959e16d97SJunchao Zhang         }                                                                        \
15059e16d97SJunchao Zhang       }                                                                          \
15159e16d97SJunchao Zhang     } else {                                                                     \
15259e16d97SJunchao Zhang       p     = MEDIAN(X,hi);                                                      \
15359e16d97SJunchao Zhang       pivot = X[p];                                                              \
15459e16d97SJunchao Zhang       TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,0,hi,l,r);                           \
15559e16d97SJunchao Zhang       ierr  = FuncName(l,X,Y,Z);CHKERRQ(ierr);                                   \
15659e16d97SJunchao Zhang       ierr  = FuncName(hi-r+1,X+r,Y+r,Z+r);CHKERRQ(ierr);                        \
15759e16d97SJunchao Zhang     }                                                                            \
15859e16d97SJunchao Zhang   } while(0)
159e5c89e4eSSatish Balay 
160e5c89e4eSSatish Balay /*@
161e5c89e4eSSatish Balay    PetscSortInt - Sorts an array of integers in place in increasing order.
162e5c89e4eSSatish Balay 
163e5c89e4eSSatish Balay    Not Collective
164e5c89e4eSSatish Balay 
165e5c89e4eSSatish Balay    Input Parameters:
166e5c89e4eSSatish Balay +  n  - number of values
16759e16d97SJunchao Zhang -  X  - array of integers
168e5c89e4eSSatish Balay 
169e5c89e4eSSatish Balay    Level: intermediate
170e5c89e4eSSatish Balay 
171e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntWithPermutation()
172e5c89e4eSSatish Balay @*/
173b4f531b9SJunchao Zhang PetscErrorCode  PetscSortInt(PetscInt n,PetscInt X[])
174e5c89e4eSSatish Balay {
17559e16d97SJunchao Zhang   PetscErrorCode ierr;
17659e16d97SJunchao Zhang   PetscInt       pivot,t1;
177e5c89e4eSSatish Balay 
178e5c89e4eSSatish Balay   PetscFunctionBegin;
17959e16d97SJunchao Zhang   QuickSort1(PetscSortInt,X,n,pivot,t1,ierr);
180e5c89e4eSSatish Balay   PetscFunctionReturn(0);
181e5c89e4eSSatish Balay }
182e5c89e4eSSatish Balay 
1831db36a52SBarry Smith /*@
18422ab5688SLisandro Dalcin    PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted input array
18522ab5688SLisandro Dalcin 
18622ab5688SLisandro Dalcin    Not Collective
18722ab5688SLisandro Dalcin 
18822ab5688SLisandro Dalcin    Input Parameters:
18922ab5688SLisandro Dalcin +  n  - number of values
19059e16d97SJunchao Zhang -  X  - sorted array of integers
19122ab5688SLisandro Dalcin 
19222ab5688SLisandro Dalcin    Output Parameter:
19322ab5688SLisandro Dalcin .  n - number of non-redundant values
19422ab5688SLisandro Dalcin 
19522ab5688SLisandro Dalcin    Level: intermediate
19622ab5688SLisandro Dalcin 
19722ab5688SLisandro Dalcin .seealso: PetscSortInt()
19822ab5688SLisandro Dalcin @*/
19959e16d97SJunchao Zhang PetscErrorCode  PetscSortedRemoveDupsInt(PetscInt *n,PetscInt X[])
20022ab5688SLisandro Dalcin {
20122ab5688SLisandro Dalcin   PetscInt i,s = 0,N = *n, b = 0;
20222ab5688SLisandro Dalcin 
20322ab5688SLisandro Dalcin   PetscFunctionBegin;
20422ab5688SLisandro Dalcin   for (i=0; i<N-1; i++) {
20559e16d97SJunchao Zhang     if (PetscUnlikely(X[b+s+1] < X[b])) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Input array is not sorted");
20659e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
20759e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
20822ab5688SLisandro Dalcin     } else s++;
20922ab5688SLisandro Dalcin   }
21022ab5688SLisandro Dalcin   *n = N - s;
21122ab5688SLisandro Dalcin   PetscFunctionReturn(0);
21222ab5688SLisandro Dalcin }
21322ab5688SLisandro Dalcin 
21422ab5688SLisandro Dalcin /*@
2151db36a52SBarry Smith    PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries
2161db36a52SBarry Smith 
2171db36a52SBarry Smith    Not Collective
2181db36a52SBarry Smith 
2191db36a52SBarry Smith    Input Parameters:
2201db36a52SBarry Smith +  n  - number of values
22159e16d97SJunchao Zhang -  X  - array of integers
2221db36a52SBarry Smith 
2231db36a52SBarry Smith    Output Parameter:
2241db36a52SBarry Smith .  n - number of non-redundant values
2251db36a52SBarry Smith 
2261db36a52SBarry Smith    Level: intermediate
2271db36a52SBarry Smith 
22822ab5688SLisandro Dalcin .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortedRemoveDupsInt()
2291db36a52SBarry Smith @*/
23059e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsInt(PetscInt *n,PetscInt X[])
2311db36a52SBarry Smith {
2321db36a52SBarry Smith   PetscErrorCode ierr;
2331db36a52SBarry Smith 
2341db36a52SBarry Smith   PetscFunctionBegin;
23559e16d97SJunchao Zhang   ierr = PetscSortInt(*n,X);CHKERRQ(ierr);
23659e16d97SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(n,X);CHKERRQ(ierr);
2371db36a52SBarry Smith   PetscFunctionReturn(0);
2381db36a52SBarry Smith }
2391db36a52SBarry Smith 
24060e03357SMatthew G Knepley /*@
241d9f1162dSJed Brown   PetscFindInt - Finds integer in a sorted array of integers
24260e03357SMatthew G Knepley 
24360e03357SMatthew G Knepley    Not Collective
24460e03357SMatthew G Knepley 
24560e03357SMatthew G Knepley    Input Parameters:
24660e03357SMatthew G Knepley +  key - the integer to locate
247d9f1162dSJed Brown .  n   - number of values in the array
24859e16d97SJunchao Zhang -  X  - array of integers
24960e03357SMatthew G Knepley 
25060e03357SMatthew G Knepley    Output Parameter:
251d9f1162dSJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
25260e03357SMatthew G Knepley 
25360e03357SMatthew G Knepley    Level: intermediate
25460e03357SMatthew G Knepley 
25560e03357SMatthew G Knepley .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
25660e03357SMatthew G Knepley @*/
25759e16d97SJunchao Zhang PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
25860e03357SMatthew G Knepley {
259d9f1162dSJed Brown   PetscInt lo = 0,hi = n;
26060e03357SMatthew G Knepley 
26160e03357SMatthew G Knepley   PetscFunctionBegin;
26260e03357SMatthew G Knepley   PetscValidPointer(loc,4);
26398781d82SMatthew G Knepley   if (!n) {*loc = -1; PetscFunctionReturn(0);}
26459e16d97SJunchao Zhang   PetscValidPointer(X,3);
265d9f1162dSJed Brown   while (hi - lo > 1) {
266d9f1162dSJed Brown     PetscInt mid = lo + (hi - lo)/2;
26759e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
268d9f1162dSJed Brown     else               lo = mid;
26960e03357SMatthew G Knepley   }
27059e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
27160e03357SMatthew G Knepley   PetscFunctionReturn(0);
27260e03357SMatthew G Knepley }
27360e03357SMatthew G Knepley 
274d2aeb606SJed Brown /*@
275*f1cab4e1SJunchao Zhang   PetscCheckDupsInt - Checks if an integer array has duplicates
276*f1cab4e1SJunchao Zhang 
277*f1cab4e1SJunchao Zhang    Not Collective
278*f1cab4e1SJunchao Zhang 
279*f1cab4e1SJunchao Zhang    Input Parameters:
280*f1cab4e1SJunchao Zhang +  n  - number of values in the array
281*f1cab4e1SJunchao Zhang -  X  - array of integers
282*f1cab4e1SJunchao Zhang 
283*f1cab4e1SJunchao Zhang 
284*f1cab4e1SJunchao Zhang    Output Parameter:
285*f1cab4e1SJunchao Zhang .  dups - True if the array has dups, otherwise false
286*f1cab4e1SJunchao Zhang 
287*f1cab4e1SJunchao Zhang    Level: intermediate
288*f1cab4e1SJunchao Zhang 
289*f1cab4e1SJunchao Zhang .seealso: PetscSortRemoveDupsInt()
290*f1cab4e1SJunchao Zhang @*/
291*f1cab4e1SJunchao Zhang PetscErrorCode PetscCheckDupsInt(PetscInt n,const PetscInt X[],PetscBool *dups)
292*f1cab4e1SJunchao Zhang {
293*f1cab4e1SJunchao Zhang   PetscErrorCode ierr;
294*f1cab4e1SJunchao Zhang   PetscInt       i;
295*f1cab4e1SJunchao Zhang   PetscHSetI     ht;
296*f1cab4e1SJunchao Zhang   PetscBool      missing;
297*f1cab4e1SJunchao Zhang 
298*f1cab4e1SJunchao Zhang   PetscFunctionBegin;
299*f1cab4e1SJunchao Zhang   PetscValidPointer(dups,3);
300*f1cab4e1SJunchao Zhang   *dups = PETSC_FALSE;
301*f1cab4e1SJunchao Zhang   if (n > 1) {
302*f1cab4e1SJunchao Zhang     ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
303*f1cab4e1SJunchao Zhang     ierr = PetscHSetIResize(ht,n);CHKERRQ(ierr);
304*f1cab4e1SJunchao Zhang     for (i=0; i<n; i++) {
305*f1cab4e1SJunchao Zhang       ierr = PetscHSetIQueryAdd(ht,X[i],&missing);CHKERRQ(ierr);
306*f1cab4e1SJunchao Zhang       if(!missing) {*dups = PETSC_TRUE; break;}
307*f1cab4e1SJunchao Zhang     }
308*f1cab4e1SJunchao Zhang     ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
309*f1cab4e1SJunchao Zhang   }
310*f1cab4e1SJunchao Zhang   PetscFunctionReturn(0);
311*f1cab4e1SJunchao Zhang }
312*f1cab4e1SJunchao Zhang 
313*f1cab4e1SJunchao Zhang /*@
314d2aeb606SJed Brown   PetscFindMPIInt - Finds MPI integer in a sorted array of integers
315d2aeb606SJed Brown 
316d2aeb606SJed Brown    Not Collective
317d2aeb606SJed Brown 
318d2aeb606SJed Brown    Input Parameters:
319d2aeb606SJed Brown +  key - the integer to locate
320d2aeb606SJed Brown .  n   - number of values in the array
32159e16d97SJunchao Zhang -  X   - array of integers
322d2aeb606SJed Brown 
323d2aeb606SJed Brown    Output Parameter:
324d2aeb606SJed Brown .  loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
325d2aeb606SJed Brown 
326d2aeb606SJed Brown    Level: intermediate
327d2aeb606SJed Brown 
328d2aeb606SJed Brown .seealso: PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
329d2aeb606SJed Brown @*/
33059e16d97SJunchao Zhang PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
331d2aeb606SJed Brown {
332d2aeb606SJed Brown   PetscInt lo = 0,hi = n;
333d2aeb606SJed Brown 
334d2aeb606SJed Brown   PetscFunctionBegin;
335d2aeb606SJed Brown   PetscValidPointer(loc,4);
336d2aeb606SJed Brown   if (!n) {*loc = -1; PetscFunctionReturn(0);}
33759e16d97SJunchao Zhang   PetscValidPointer(X,3);
338d2aeb606SJed Brown   while (hi - lo > 1) {
339d2aeb606SJed Brown     PetscInt mid = lo + (hi - lo)/2;
34059e16d97SJunchao Zhang     if (key < X[mid]) hi = mid;
341d2aeb606SJed Brown     else               lo = mid;
342d2aeb606SJed Brown   }
34359e16d97SJunchao Zhang   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
344e5c89e4eSSatish Balay   PetscFunctionReturn(0);
345e5c89e4eSSatish Balay }
346e5c89e4eSSatish Balay 
347e5c89e4eSSatish Balay /*@
348e5c89e4eSSatish Balay    PetscSortIntWithArray - Sorts an array of integers in place in increasing order;
349e5c89e4eSSatish Balay        changes a second array to match the sorted first array.
350e5c89e4eSSatish Balay 
351e5c89e4eSSatish Balay    Not Collective
352e5c89e4eSSatish Balay 
353e5c89e4eSSatish Balay    Input Parameters:
354e5c89e4eSSatish Balay +  n  - number of values
35559e16d97SJunchao Zhang .  X  - array of integers
35659e16d97SJunchao Zhang -  Y  - second array of integers
357e5c89e4eSSatish Balay 
358e5c89e4eSSatish Balay    Level: intermediate
359e5c89e4eSSatish Balay 
360e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
361e5c89e4eSSatish Balay @*/
362b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithArray(PetscInt n,PetscInt X[],PetscInt Y[])
363e5c89e4eSSatish Balay {
364e5c89e4eSSatish Balay   PetscErrorCode ierr;
36559e16d97SJunchao Zhang   PetscInt       pivot,t1,t2;
366e5c89e4eSSatish Balay 
367e5c89e4eSSatish Balay   PetscFunctionBegin;
36859e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithArray,X,Y,n,pivot,t1,t2,ierr);
369c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
370c1f0200aSDmitry Karpeev }
371c1f0200aSDmitry Karpeev 
372c1f0200aSDmitry Karpeev /*@
373c1f0200aSDmitry Karpeev    PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order;
374c1f0200aSDmitry Karpeev        changes a pair of integer arrays to match the sorted first array.
375c1f0200aSDmitry Karpeev 
376c1f0200aSDmitry Karpeev    Not Collective
377c1f0200aSDmitry Karpeev 
378c1f0200aSDmitry Karpeev    Input Parameters:
379c1f0200aSDmitry Karpeev +  n  - number of values
38059e16d97SJunchao Zhang .  X  - array of integers
38159e16d97SJunchao Zhang .  Y  - second array of integers (first array of the pair)
38259e16d97SJunchao Zhang -  Z  - third array of integers  (second array of the pair)
383c1f0200aSDmitry Karpeev 
384c1f0200aSDmitry Karpeev    Level: intermediate
385c1f0200aSDmitry Karpeev 
386c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray()
387c1f0200aSDmitry Karpeev @*/
388b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithArrayPair(PetscInt n,PetscInt X[],PetscInt Y[],PetscInt Z[])
389c1f0200aSDmitry Karpeev {
390c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
39159e16d97SJunchao Zhang   PetscInt       pivot,t1,t2,t3;
392c1f0200aSDmitry Karpeev 
393c1f0200aSDmitry Karpeev   PetscFunctionBegin;
39459e16d97SJunchao Zhang   QuickSort3(PetscSortIntWithArrayPair,X,Y,Z,n,pivot,t1,t2,t3,ierr);
395c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
396c1f0200aSDmitry Karpeev }
397c1f0200aSDmitry Karpeev 
39817d7d925SStefano Zampini /*@
39917d7d925SStefano Zampini    PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order.
40017d7d925SStefano Zampini 
40117d7d925SStefano Zampini    Not Collective
40217d7d925SStefano Zampini 
40317d7d925SStefano Zampini    Input Parameters:
40417d7d925SStefano Zampini +  n  - number of values
40559e16d97SJunchao Zhang -  X  - array of integers
40617d7d925SStefano Zampini 
40717d7d925SStefano Zampini    Level: intermediate
40817d7d925SStefano Zampini 
40917d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation()
41017d7d925SStefano Zampini @*/
411b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIInt(PetscInt n,PetscMPIInt X[])
41217d7d925SStefano Zampini {
41359e16d97SJunchao Zhang   PetscErrorCode ierr;
41459e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
41517d7d925SStefano Zampini 
41617d7d925SStefano Zampini   PetscFunctionBegin;
41759e16d97SJunchao Zhang   QuickSort1(PetscSortMPIInt,X,n,pivot,t1,ierr);
41817d7d925SStefano Zampini   PetscFunctionReturn(0);
41917d7d925SStefano Zampini }
42017d7d925SStefano Zampini 
42117d7d925SStefano Zampini /*@
42217d7d925SStefano Zampini    PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries
42317d7d925SStefano Zampini 
42417d7d925SStefano Zampini    Not Collective
42517d7d925SStefano Zampini 
42617d7d925SStefano Zampini    Input Parameters:
42717d7d925SStefano Zampini +  n  - number of values
42859e16d97SJunchao Zhang -  X  - array of integers
42917d7d925SStefano Zampini 
43017d7d925SStefano Zampini    Output Parameter:
43117d7d925SStefano Zampini .  n - number of non-redundant values
43217d7d925SStefano Zampini 
43317d7d925SStefano Zampini    Level: intermediate
43417d7d925SStefano Zampini 
43517d7d925SStefano Zampini .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
43617d7d925SStefano Zampini @*/
43759e16d97SJunchao Zhang PetscErrorCode  PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt X[])
43817d7d925SStefano Zampini {
43917d7d925SStefano Zampini   PetscErrorCode ierr;
44017d7d925SStefano Zampini   PetscInt       i,s = 0,N = *n, b = 0;
44117d7d925SStefano Zampini 
44217d7d925SStefano Zampini   PetscFunctionBegin;
44359e16d97SJunchao Zhang   ierr = PetscSortMPIInt(N,X);CHKERRQ(ierr);
44417d7d925SStefano Zampini   for (i=0; i<N-1; i++) {
44559e16d97SJunchao Zhang     if (X[b+s+1] != X[b]) {
44659e16d97SJunchao Zhang       X[b+1] = X[b+s+1]; b++;
447a297a907SKarl Rupp     } else s++;
44817d7d925SStefano Zampini   }
44917d7d925SStefano Zampini   *n = N - s;
45017d7d925SStefano Zampini   PetscFunctionReturn(0);
45117d7d925SStefano Zampini }
452c1f0200aSDmitry Karpeev 
4534d615ea0SBarry Smith /*@
4544d615ea0SBarry Smith    PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order;
4554d615ea0SBarry Smith        changes a second array to match the sorted first array.
4564d615ea0SBarry Smith 
4574d615ea0SBarry Smith    Not Collective
4584d615ea0SBarry Smith 
4594d615ea0SBarry Smith    Input Parameters:
4604d615ea0SBarry Smith +  n  - number of values
46159e16d97SJunchao Zhang .  X  - array of integers
46259e16d97SJunchao Zhang -  Y  - second array of integers
4634d615ea0SBarry Smith 
4644d615ea0SBarry Smith    Level: intermediate
4654d615ea0SBarry Smith 
4664d615ea0SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
4674d615ea0SBarry Smith @*/
468b4f531b9SJunchao Zhang PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt X[],PetscMPIInt Y[])
4694d615ea0SBarry Smith {
4704d615ea0SBarry Smith   PetscErrorCode ierr;
47159e16d97SJunchao Zhang   PetscMPIInt    pivot,t1,t2;
4724d615ea0SBarry Smith 
4734d615ea0SBarry Smith   PetscFunctionBegin;
47459e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithArray,X,Y,n,pivot,t1,t2,ierr);
4755569a785SJunchao Zhang   PetscFunctionReturn(0);
4765569a785SJunchao Zhang }
4775569a785SJunchao Zhang 
4785569a785SJunchao Zhang /*@
4795569a785SJunchao Zhang    PetscSortMPIIntWithIntArray - Sorts an array of MPI integers in place in increasing order;
4805569a785SJunchao Zhang        changes a second array of Petsc intergers to match the sorted first array.
4815569a785SJunchao Zhang 
4825569a785SJunchao Zhang    Not Collective
4835569a785SJunchao Zhang 
4845569a785SJunchao Zhang    Input Parameters:
4855569a785SJunchao Zhang +  n  - number of values
48659e16d97SJunchao Zhang .  X  - array of MPI integers
48759e16d97SJunchao Zhang -  Y  - second array of Petsc integers
4885569a785SJunchao Zhang 
4895569a785SJunchao Zhang    Level: intermediate
4905569a785SJunchao Zhang 
4915569a785SJunchao Zhang    Notes: this routine is useful when one needs to sort MPI ranks with other integer arrays.
4925569a785SJunchao Zhang 
4935569a785SJunchao Zhang .seealso: PetscSortMPIIntWithArray()
4945569a785SJunchao Zhang @*/
495b4f531b9SJunchao Zhang PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n,PetscMPIInt X[],PetscInt Y[])
4965569a785SJunchao Zhang {
4975569a785SJunchao Zhang   PetscErrorCode ierr;
49859e16d97SJunchao Zhang   PetscMPIInt    pivot,t1;
4995569a785SJunchao Zhang   PetscInt       t2;
5005569a785SJunchao Zhang 
5015569a785SJunchao Zhang   PetscFunctionBegin;
50259e16d97SJunchao Zhang   QuickSort2(PetscSortMPIIntWithIntArray,X,Y,n,pivot,t1,t2,ierr);
503e5c89e4eSSatish Balay   PetscFunctionReturn(0);
504e5c89e4eSSatish Balay }
505e5c89e4eSSatish Balay 
506e5c89e4eSSatish Balay /*@
507e5c89e4eSSatish Balay    PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
508e5c89e4eSSatish Balay        changes a second SCALAR array to match the sorted first INTEGER array.
509e5c89e4eSSatish Balay 
510e5c89e4eSSatish Balay    Not Collective
511e5c89e4eSSatish Balay 
512e5c89e4eSSatish Balay    Input Parameters:
513e5c89e4eSSatish Balay +  n  - number of values
51459e16d97SJunchao Zhang .  X  - array of integers
51559e16d97SJunchao Zhang -  Y  - second array of scalars
516e5c89e4eSSatish Balay 
517e5c89e4eSSatish Balay    Level: intermediate
518e5c89e4eSSatish Balay 
519e5c89e4eSSatish Balay .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
520e5c89e4eSSatish Balay @*/
521b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithScalarArray(PetscInt n,PetscInt X[],PetscScalar Y[])
522e5c89e4eSSatish Balay {
523e5c89e4eSSatish Balay   PetscErrorCode ierr;
52459e16d97SJunchao Zhang   PetscInt       pivot,t1;
52559e16d97SJunchao Zhang   PetscScalar    t2;
526e5c89e4eSSatish Balay 
527e5c89e4eSSatish Balay   PetscFunctionBegin;
52859e16d97SJunchao Zhang   QuickSort2(PetscSortIntWithScalarArray,X,Y,n,pivot,t1,t2,ierr);
52917df18f2SToby Isaac   PetscFunctionReturn(0);
53017df18f2SToby Isaac }
53117df18f2SToby Isaac 
5326c2863d0SBarry Smith /*@C
53317df18f2SToby Isaac    PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order;
53417df18f2SToby Isaac        changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
53517df18f2SToby Isaac        provide workspace (the size of an element in the data array) to use when sorting.
53617df18f2SToby Isaac 
53717df18f2SToby Isaac    Not Collective
53817df18f2SToby Isaac 
53917df18f2SToby Isaac    Input Parameters:
54017df18f2SToby Isaac +  n  - number of values
54159e16d97SJunchao Zhang .  X  - array of integers
54259e16d97SJunchao Zhang .  Y  - second array of data
54317df18f2SToby Isaac .  size - sizeof elements in the data array in bytes
54459e16d97SJunchao Zhang -  t2   - workspace of "size" bytes used when sorting
54517df18f2SToby Isaac 
54617df18f2SToby Isaac    Level: intermediate
54717df18f2SToby Isaac 
54817df18f2SToby Isaac .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
54917df18f2SToby Isaac @*/
550b4f531b9SJunchao Zhang PetscErrorCode  PetscSortIntWithDataArray(PetscInt n,PetscInt X[],void *Y,size_t size,void *t2)
55117df18f2SToby Isaac {
55217df18f2SToby Isaac   PetscErrorCode ierr;
55359e16d97SJunchao Zhang   char           *YY = (char*)Y;
55459e16d97SJunchao Zhang   PetscInt       i,j,p,t1,pivot,hi=n-1,l,r;
55517df18f2SToby Isaac 
55617df18f2SToby Isaac   PetscFunctionBegin;
55717df18f2SToby Isaac   if (n<8) {
55859e16d97SJunchao Zhang     for (i=0; i<n; i++) {
55959e16d97SJunchao Zhang       pivot = X[i];
56059e16d97SJunchao Zhang       for (j=i+1; j<n; j++) {
56159e16d97SJunchao Zhang         if (pivot > X[j]) {
56259e16d97SJunchao Zhang           SWAP2Data(X[i],X[j],YY+size*i,YY+size*j,t1,t2,size);
56359e16d97SJunchao Zhang           pivot = X[i];
56417df18f2SToby Isaac         }
56517df18f2SToby Isaac       }
56617df18f2SToby Isaac     }
56717df18f2SToby Isaac   } else {
56859e16d97SJunchao Zhang     /* Two way partition */
56959e16d97SJunchao Zhang     p     = MEDIAN(X,hi);
57059e16d97SJunchao Zhang     pivot = X[p];
57159e16d97SJunchao Zhang     l     = 0;
57259e16d97SJunchao Zhang     r     = hi;
57359e16d97SJunchao Zhang     while(1) {
57459e16d97SJunchao Zhang       while (X[l] < pivot) l++;
57559e16d97SJunchao Zhang       while (X[r] > pivot) r--;
57659e16d97SJunchao Zhang       if (l >= r) {r++; break;}
57759e16d97SJunchao Zhang       SWAP2Data(X[l],X[r],YY+size*l,YY+size*r,t1,t2,size);
57859e16d97SJunchao Zhang       l++;
57959e16d97SJunchao Zhang       r--;
58059e16d97SJunchao Zhang     }
58159e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(l,X,Y,size,t2);CHKERRQ(ierr);
58259e16d97SJunchao Zhang     ierr = PetscSortIntWithDataArray(hi-r+1,X+r,YY+size*r,size,t2);CHKERRQ(ierr);
58317df18f2SToby Isaac   }
58417df18f2SToby Isaac   PetscFunctionReturn(0);
58517df18f2SToby Isaac }
58617df18f2SToby Isaac 
58721e72a00SBarry Smith /*@
58821e72a00SBarry Smith    PetscMergeIntArray -     Merges two SORTED integer arrays, removes duplicate elements.
58921e72a00SBarry Smith 
59021e72a00SBarry Smith    Not Collective
59121e72a00SBarry Smith 
59221e72a00SBarry Smith    Input Parameters:
59321e72a00SBarry Smith +  an  - number of values in the first array
59421e72a00SBarry Smith .  aI  - first sorted array of integers
59521e72a00SBarry Smith .  bn  - number of values in the second array
59621e72a00SBarry Smith -  bI  - second array of integers
59721e72a00SBarry Smith 
59821e72a00SBarry Smith    Output Parameters:
59921e72a00SBarry Smith +  n   - number of values in the merged array
6006c2863d0SBarry Smith -  L   - merged sorted array, this is allocated if an array is not provided
60121e72a00SBarry Smith 
60221e72a00SBarry Smith    Level: intermediate
60321e72a00SBarry Smith 
60421e72a00SBarry Smith .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
60521e72a00SBarry Smith @*/
6066c2863d0SBarry Smith PetscErrorCode  PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
60721e72a00SBarry Smith {
60821e72a00SBarry Smith   PetscErrorCode ierr;
60921e72a00SBarry Smith   PetscInt       *L_ = *L, ak, bk, k;
61021e72a00SBarry Smith 
61121e72a00SBarry Smith   if (!L_) {
61221e72a00SBarry Smith     ierr = PetscMalloc1(an+bn, L);CHKERRQ(ierr);
61321e72a00SBarry Smith     L_   = *L;
61421e72a00SBarry Smith   }
61521e72a00SBarry Smith   k = ak = bk = 0;
61621e72a00SBarry Smith   while (ak < an && bk < bn) {
61721e72a00SBarry Smith     if (aI[ak] == bI[bk]) {
61821e72a00SBarry Smith       L_[k] = aI[ak];
61921e72a00SBarry Smith       ++ak;
62021e72a00SBarry Smith       ++bk;
62121e72a00SBarry Smith       ++k;
62221e72a00SBarry Smith     } else if (aI[ak] < bI[bk]) {
62321e72a00SBarry Smith       L_[k] = aI[ak];
62421e72a00SBarry Smith       ++ak;
62521e72a00SBarry Smith       ++k;
62621e72a00SBarry Smith     } else {
62721e72a00SBarry Smith       L_[k] = bI[bk];
62821e72a00SBarry Smith       ++bk;
62921e72a00SBarry Smith       ++k;
63021e72a00SBarry Smith     }
63121e72a00SBarry Smith   }
63221e72a00SBarry Smith   if (ak < an) {
633580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
63421e72a00SBarry Smith     k   += (an-ak);
63521e72a00SBarry Smith   }
63621e72a00SBarry Smith   if (bk < bn) {
637580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
63821e72a00SBarry Smith     k   += (bn-bk);
63921e72a00SBarry Smith   }
64021e72a00SBarry Smith   *n = k;
64121e72a00SBarry Smith   PetscFunctionReturn(0);
64221e72a00SBarry Smith }
643b4301105SBarry Smith 
644c1f0200aSDmitry Karpeev /*@
64521e72a00SBarry Smith    PetscMergeIntArrayPair -     Merges two SORTED integer arrays that share NO common values along with an additional array of integers.
646c1f0200aSDmitry Karpeev                                 The additional arrays are the same length as sorted arrays and are merged
647c1f0200aSDmitry Karpeev                                 in the order determined by the merging of the sorted pair.
648c1f0200aSDmitry Karpeev 
649c1f0200aSDmitry Karpeev    Not Collective
650c1f0200aSDmitry Karpeev 
651c1f0200aSDmitry Karpeev    Input Parameters:
652c1f0200aSDmitry Karpeev +  an  - number of values in the first array
653c1f0200aSDmitry Karpeev .  aI  - first sorted array of integers
654c1f0200aSDmitry Karpeev .  aJ  - first additional array of integers
655c1f0200aSDmitry Karpeev .  bn  - number of values in the second array
656c1f0200aSDmitry Karpeev .  bI  - second array of integers
657c1f0200aSDmitry Karpeev -  bJ  - second additional of integers
658c1f0200aSDmitry Karpeev 
659c1f0200aSDmitry Karpeev    Output Parameters:
660c1f0200aSDmitry Karpeev +  n   - number of values in the merged array (== an + bn)
66114c49006SJed Brown .  L   - merged sorted array
662c1f0200aSDmitry Karpeev -  J   - merged additional array
663c1f0200aSDmitry Karpeev 
66495452b02SPatrick Sanan    Notes:
66595452b02SPatrick 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
666c1f0200aSDmitry Karpeev    Level: intermediate
667c1f0200aSDmitry Karpeev 
668c1f0200aSDmitry Karpeev .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
669c1f0200aSDmitry Karpeev @*/
6706c2863d0SBarry 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)
671c1f0200aSDmitry Karpeev {
672c1f0200aSDmitry Karpeev   PetscErrorCode ierr;
67370d8d27cSBarry Smith   PetscInt       n_, *L_, *J_, ak, bk, k;
674c1f0200aSDmitry Karpeev 
67514c49006SJed Brown   PetscFunctionBegin;
67670d8d27cSBarry Smith   PetscValidIntPointer(L,8);
67770d8d27cSBarry Smith   PetscValidIntPointer(J,9);
678c1f0200aSDmitry Karpeev   n_ = an + bn;
679c1f0200aSDmitry Karpeev   *n = n_;
68070d8d27cSBarry Smith   if (!*L) {
681785e854fSJed Brown     ierr = PetscMalloc1(n_, L);CHKERRQ(ierr);
68270d8d27cSBarry Smith   }
683d7aa01a8SHong Zhang   L_ = *L;
68470d8d27cSBarry Smith   if (!*J) {
68570d8d27cSBarry Smith     ierr = PetscMalloc1(n_, J);CHKERRQ(ierr);
686c1f0200aSDmitry Karpeev   }
687c1f0200aSDmitry Karpeev   J_   = *J;
688c1f0200aSDmitry Karpeev   k = ak = bk = 0;
689c1f0200aSDmitry Karpeev   while (ak < an && bk < bn) {
690c1f0200aSDmitry Karpeev     if (aI[ak] <= bI[bk]) {
691d7aa01a8SHong Zhang       L_[k] = aI[ak];
692c1f0200aSDmitry Karpeev       J_[k] = aJ[ak];
693c1f0200aSDmitry Karpeev       ++ak;
694c1f0200aSDmitry Karpeev       ++k;
695a297a907SKarl Rupp     } else {
696d7aa01a8SHong Zhang       L_[k] = bI[bk];
697c1f0200aSDmitry Karpeev       J_[k] = bJ[bk];
698c1f0200aSDmitry Karpeev       ++bk;
699c1f0200aSDmitry Karpeev       ++k;
700c1f0200aSDmitry Karpeev     }
701c1f0200aSDmitry Karpeev   }
702c1f0200aSDmitry Karpeev   if (ak < an) {
703580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,aI+ak,an-ak);CHKERRQ(ierr);
704580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,aJ+ak,an-ak);CHKERRQ(ierr);
705c1f0200aSDmitry Karpeev     k   += (an-ak);
706c1f0200aSDmitry Karpeev   }
707c1f0200aSDmitry Karpeev   if (bk < bn) {
708580bdb30SBarry Smith     ierr = PetscArraycpy(L_+k,bI+bk,bn-bk);CHKERRQ(ierr);
709580bdb30SBarry Smith     ierr = PetscArraycpy(J_+k,bJ+bk,bn-bk);CHKERRQ(ierr);
710c1f0200aSDmitry Karpeev   }
711c1f0200aSDmitry Karpeev   PetscFunctionReturn(0);
712c1f0200aSDmitry Karpeev }
713c1f0200aSDmitry Karpeev 
714e498c390SJed Brown /*@
715e498c390SJed Brown    PetscMergeMPIIntArray -     Merges two SORTED integer arrays.
716e498c390SJed Brown 
717e498c390SJed Brown    Not Collective
718e498c390SJed Brown 
719e498c390SJed Brown    Input Parameters:
720e498c390SJed Brown +  an  - number of values in the first array
721e498c390SJed Brown .  aI  - first sorted array of integers
722e498c390SJed Brown .  bn  - number of values in the second array
723e498c390SJed Brown -  bI  - second array of integers
724e498c390SJed Brown 
725e498c390SJed Brown    Output Parameters:
726e498c390SJed Brown +  n   - number of values in the merged array (<= an + bn)
727e498c390SJed Brown -  L   - merged sorted array, allocated if address of NULL pointer is passed
728e498c390SJed Brown 
729e498c390SJed Brown    Level: intermediate
730e498c390SJed Brown 
731e498c390SJed Brown .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
732e498c390SJed Brown @*/
733e498c390SJed Brown PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L)
734e498c390SJed Brown {
735e498c390SJed Brown   PetscErrorCode ierr;
736e498c390SJed Brown   PetscInt       ai,bi,k;
737e498c390SJed Brown 
738e498c390SJed Brown   PetscFunctionBegin;
739e498c390SJed Brown   if (!*L) {ierr = PetscMalloc1((an+bn),L);CHKERRQ(ierr);}
740e498c390SJed Brown   for (ai=0,bi=0,k=0; ai<an || bi<bn; ) {
741e498c390SJed Brown     PetscInt t = -1;
742e498c390SJed Brown     for ( ; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
743e498c390SJed Brown     for ( ; bi<bn && bI[bi] == t; bi++);
744e498c390SJed Brown     for ( ; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
745e498c390SJed Brown     for ( ; ai<an && aI[ai] == t; ai++);
746e498c390SJed Brown   }
747e498c390SJed Brown   *n = k;
748e498c390SJed Brown   PetscFunctionReturn(0);
749e498c390SJed Brown }
750e5c89e4eSSatish Balay 
7516c2863d0SBarry Smith /*@C
75242eaa7f3SBarry Smith    PetscProcessTree - Prepares tree data to be displayed graphically
75342eaa7f3SBarry Smith 
75442eaa7f3SBarry Smith    Not Collective
75542eaa7f3SBarry Smith 
75642eaa7f3SBarry Smith    Input Parameters:
75742eaa7f3SBarry Smith +  n  - number of values
75842eaa7f3SBarry Smith .  mask - indicates those entries in the tree, location 0 is always masked
75942eaa7f3SBarry Smith -  parentid - indicates the parent of each entry
76042eaa7f3SBarry Smith 
76142eaa7f3SBarry Smith    Output Parameters:
76242eaa7f3SBarry Smith +  Nlevels - the number of levels
76342eaa7f3SBarry Smith .  Level - for each node tells its level
76442eaa7f3SBarry Smith .  Levelcnts - the number of nodes on each level
76542eaa7f3SBarry Smith .  Idbylevel - a list of ids on each of the levels, first level followed by second etc
76642eaa7f3SBarry Smith -  Column - for each id tells its column index
76742eaa7f3SBarry Smith 
7686c2863d0SBarry Smith    Level: developer
76942eaa7f3SBarry Smith 
77095452b02SPatrick Sanan    Notes:
77195452b02SPatrick Sanan     This code is not currently used
77242eaa7f3SBarry Smith 
77342eaa7f3SBarry Smith .seealso: PetscSortReal(), PetscSortIntWithPermutation()
77442eaa7f3SBarry Smith @*/
7757087cfbeSBarry Smith PetscErrorCode  PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column)
77642eaa7f3SBarry Smith {
77742eaa7f3SBarry Smith   PetscInt       i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column;
77842eaa7f3SBarry Smith   PetscErrorCode ierr;
779ace3abfcSBarry Smith   PetscBool      done = PETSC_FALSE;
78042eaa7f3SBarry Smith 
78142eaa7f3SBarry Smith   PetscFunctionBegin;
78242eaa7f3SBarry Smith   if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set");
78342eaa7f3SBarry Smith   for (i=0; i<n; i++) {
78442eaa7f3SBarry Smith     if (mask[i]) continue;
78542eaa7f3SBarry Smith     if (parentid[i]  == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent");
78642eaa7f3SBarry Smith     if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked");
78742eaa7f3SBarry Smith   }
78842eaa7f3SBarry Smith 
78942eaa7f3SBarry Smith   for (i=0; i<n; i++) {
79042eaa7f3SBarry Smith     if (!mask[i]) nmask++;
79142eaa7f3SBarry Smith   }
79242eaa7f3SBarry Smith 
79342eaa7f3SBarry Smith   /* determine the level in the tree of each node */
7941795a4d1SJed Brown   ierr = PetscCalloc1(n,&level);CHKERRQ(ierr);
795a297a907SKarl Rupp 
79642eaa7f3SBarry Smith   level[0] = 1;
79742eaa7f3SBarry Smith   while (!done) {
79842eaa7f3SBarry Smith     done = PETSC_TRUE;
79942eaa7f3SBarry Smith     for (i=0; i<n; i++) {
80042eaa7f3SBarry Smith       if (mask[i]) continue;
80142eaa7f3SBarry Smith       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
80242eaa7f3SBarry Smith       else if (!level[i]) done = PETSC_FALSE;
80342eaa7f3SBarry Smith     }
80442eaa7f3SBarry Smith   }
80542eaa7f3SBarry Smith   for (i=0; i<n; i++) {
80642eaa7f3SBarry Smith     level[i]--;
80742eaa7f3SBarry Smith     nlevels = PetscMax(nlevels,level[i]);
80842eaa7f3SBarry Smith   }
80942eaa7f3SBarry Smith 
81042eaa7f3SBarry Smith   /* count the number of nodes on each level and its max */
8111795a4d1SJed Brown   ierr = PetscCalloc1(nlevels,&levelcnt);CHKERRQ(ierr);
81242eaa7f3SBarry Smith   for (i=0; i<n; i++) {
81342eaa7f3SBarry Smith     if (mask[i]) continue;
81442eaa7f3SBarry Smith     levelcnt[level[i]-1]++;
81542eaa7f3SBarry Smith   }
816a297a907SKarl Rupp   for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]);
81742eaa7f3SBarry Smith 
81842eaa7f3SBarry Smith   /* for each level sort the ids by the parent id */
819dcca6d9dSJed Brown   ierr = PetscMalloc2(levelmax,&workid,levelmax,&workparentid);CHKERRQ(ierr);
820785e854fSJed Brown   ierr = PetscMalloc1(nmask,&idbylevel);CHKERRQ(ierr);
82142eaa7f3SBarry Smith   for (j=1; j<=nlevels;j++) {
82242eaa7f3SBarry Smith     cnt = 0;
82342eaa7f3SBarry Smith     for (i=0; i<n; i++) {
82442eaa7f3SBarry Smith       if (mask[i]) continue;
82542eaa7f3SBarry Smith       if (level[i] != j) continue;
82642eaa7f3SBarry Smith       workid[cnt]         = i;
82742eaa7f3SBarry Smith       workparentid[cnt++] = parentid[i];
82842eaa7f3SBarry Smith     }
82942eaa7f3SBarry Smith     /*  PetscIntView(cnt,workparentid,0);
83042eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);
83142eaa7f3SBarry Smith     ierr = PetscSortIntWithArray(cnt,workparentid,workid);CHKERRQ(ierr);
83242eaa7f3SBarry Smith     PetscIntView(cnt,workparentid,0);
83342eaa7f3SBarry Smith     PetscIntView(cnt,workid,0);*/
834580bdb30SBarry Smith     ierr  = PetscArraycpy(idbylevel+tcnt,workid,cnt);CHKERRQ(ierr);
83542eaa7f3SBarry Smith     tcnt += cnt;
83642eaa7f3SBarry Smith   }
83742eaa7f3SBarry Smith   if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes");
83842eaa7f3SBarry Smith   ierr = PetscFree2(workid,workparentid);CHKERRQ(ierr);
83942eaa7f3SBarry Smith 
84042eaa7f3SBarry Smith   /* for each node list its column */
841785e854fSJed Brown   ierr = PetscMalloc1(n,&column);CHKERRQ(ierr);
84242eaa7f3SBarry Smith   cnt = 0;
84342eaa7f3SBarry Smith   for (j=0; j<nlevels; j++) {
84442eaa7f3SBarry Smith     for (i=0; i<levelcnt[j]; i++) {
84542eaa7f3SBarry Smith       column[idbylevel[cnt++]] = i;
84642eaa7f3SBarry Smith     }
84742eaa7f3SBarry Smith   }
84842eaa7f3SBarry Smith 
84942eaa7f3SBarry Smith   *Nlevels   = nlevels;
85042eaa7f3SBarry Smith   *Level     = level;
85142eaa7f3SBarry Smith   *Levelcnt  = levelcnt;
85242eaa7f3SBarry Smith   *Idbylevel = idbylevel;
85342eaa7f3SBarry Smith   *Column    = column;
85442eaa7f3SBarry Smith   PetscFunctionReturn(0);
85542eaa7f3SBarry Smith }
856