xref: /petsc/src/sys/utils/sortip.c (revision eca7e54bf87a86b0ac7ad3e6d18457430ad72719)
1e5c89e4eSSatish Balay /*
2e5c89e4eSSatish Balay    This file contains routines for sorting integers and doubles with a permutation array.
3e5c89e4eSSatish Balay  */
4*c435d973SPierre Jolivet #include <petsc/private/petscimpl.h>
5c6db04a5SJed Brown #include <petscsys.h> /*I  "petscsys.h"  I*/
6e5c89e4eSSatish Balay 
79371c9d4SSatish Balay #define SWAP(a, b, t) \
8a8f51744SPierre Jolivet   do { \
99371c9d4SSatish Balay     t = a; \
109371c9d4SSatish Balay     a = b; \
119371c9d4SSatish Balay     b = t; \
12a8f51744SPierre Jolivet   } while (0)
13e5c89e4eSSatish Balay 
14*c435d973SPierre Jolivet #if PetscDefined(USE_DEBUG)
15*c435d973SPierre Jolivet   #define PetscCheckIdentity(n, idx) \
16*c435d973SPierre Jolivet     do { \
17*c435d973SPierre Jolivet       for (PetscInt i = 0; i < n; ++i) PetscCheck(idx[i] == i, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input array needs to be initialized to 0:%" PetscInt_FMT, n - 1); \
18*c435d973SPierre Jolivet     } while (0)
19*c435d973SPierre Jolivet #else
20*c435d973SPierre Jolivet   #define PetscCheckIdentity(n, idx) (void)0
21*c435d973SPierre Jolivet #endif
22*c435d973SPierre Jolivet 
PetscSortIntWithPermutation_Private(const PetscInt v[],PetscInt vdx[],PetscInt right)23d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSortIntWithPermutation_Private(const PetscInt v[], PetscInt vdx[], PetscInt right)
24d71ae5a4SJacob Faibussowitsch {
25e5c89e4eSSatish Balay   PetscInt tmp, i, vl, last;
26e5c89e4eSSatish Balay 
27e5c89e4eSSatish Balay   PetscFunctionBegin;
28e5c89e4eSSatish Balay   if (right <= 1) {
29e5c89e4eSSatish Balay     if (right == 1) {
30e5c89e4eSSatish Balay       if (v[vdx[0]] > v[vdx[1]]) SWAP(vdx[0], vdx[1], tmp);
31e5c89e4eSSatish Balay     }
323ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
33e5c89e4eSSatish Balay   }
34e5c89e4eSSatish Balay   SWAP(vdx[0], vdx[right / 2], tmp);
35e5c89e4eSSatish Balay   vl   = v[vdx[0]];
36e5c89e4eSSatish Balay   last = 0;
37e5c89e4eSSatish Balay   for (i = 1; i <= right; i++) {
389371c9d4SSatish Balay     if (v[vdx[i]] < vl) {
399371c9d4SSatish Balay       last++;
409371c9d4SSatish Balay       SWAP(vdx[last], vdx[i], tmp);
419371c9d4SSatish Balay     }
42e5c89e4eSSatish Balay   }
43e5c89e4eSSatish Balay   SWAP(vdx[0], vdx[last], tmp);
449566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithPermutation_Private(v, vdx, last - 1));
459566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithPermutation_Private(v, vdx + last + 1, right - (last + 1)));
463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47e5c89e4eSSatish Balay }
48e5c89e4eSSatish Balay 
49e5c89e4eSSatish Balay /*@
50811af0c4SBarry Smith   PetscSortIntWithPermutation - Computes the permutation of `PetscInt` that gives
51e5c89e4eSSatish Balay   a sorted sequence.
52e5c89e4eSSatish Balay 
53e5c89e4eSSatish Balay   Not Collective
54e5c89e4eSSatish Balay 
55e5c89e4eSSatish Balay   Input Parameters:
56e5c89e4eSSatish Balay + n   - number of values to sort
57e5c89e4eSSatish Balay . i   - values to sort
58*c435d973SPierre Jolivet - idx - permutation array. Must be initialized to 0:`n`-1 on input.
59e5c89e4eSSatish Balay 
60e5c89e4eSSatish Balay   Level: intermediate
61e5c89e4eSSatish Balay 
62811af0c4SBarry Smith   Note:
63*c435d973SPierre Jolivet   On output, `i` is unchanged and `idx[j]` is the position of the `j`th smallest `PetscInt` in `i`.
64e5c89e4eSSatish Balay 
65db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()`, `PetscSortIntWithArray()`
66e5c89e4eSSatish Balay  @*/
PetscSortIntWithPermutation(PetscInt n,const PetscInt i[],PetscInt idx[])67d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortIntWithPermutation(PetscInt n, const PetscInt i[], PetscInt idx[])
68d71ae5a4SJacob Faibussowitsch {
69e5c89e4eSSatish Balay   PetscInt j, k, tmp, ik;
70e5c89e4eSSatish Balay 
71e5c89e4eSSatish Balay   PetscFunctionBegin;
72*c435d973SPierre Jolivet   if (n > 0) {
73*c435d973SPierre Jolivet     PetscAssertPointer(i, 2);
74*c435d973SPierre Jolivet     PetscAssertPointer(idx, 3);
75*c435d973SPierre Jolivet     PetscCheckIdentity(n, idx);
76*c435d973SPierre Jolivet   }
77e5c89e4eSSatish Balay   if (n < 8) {
78e5c89e4eSSatish Balay     for (k = 0; k < n; k++) {
79e5c89e4eSSatish Balay       ik = i[idx[k]];
80e5c89e4eSSatish Balay       for (j = k + 1; j < n; j++) {
81e5c89e4eSSatish Balay         if (ik > i[idx[j]]) {
82e5c89e4eSSatish Balay           SWAP(idx[k], idx[j], tmp);
83e5c89e4eSSatish Balay           ik = i[idx[k]];
84e5c89e4eSSatish Balay         }
85e5c89e4eSSatish Balay       }
86e5c89e4eSSatish Balay     }
87e5c89e4eSSatish Balay   } else {
889566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithPermutation_Private(i, idx, n - 1));
89e5c89e4eSSatish Balay   }
903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91e5c89e4eSSatish Balay }
92e5c89e4eSSatish Balay 
PetscSortRealWithPermutation_Private(const PetscReal v[],PetscInt vdx[],PetscInt right)93d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSortRealWithPermutation_Private(const PetscReal v[], PetscInt vdx[], PetscInt right)
94d71ae5a4SJacob Faibussowitsch {
95e5c89e4eSSatish Balay   PetscReal vl;
96e5c89e4eSSatish Balay   PetscInt  tmp, i, last;
97e5c89e4eSSatish Balay 
98e5c89e4eSSatish Balay   PetscFunctionBegin;
99e5c89e4eSSatish Balay   if (right <= 1) {
100e5c89e4eSSatish Balay     if (right == 1) {
101e5c89e4eSSatish Balay       if (v[vdx[0]] > v[vdx[1]]) SWAP(vdx[0], vdx[1], tmp);
102e5c89e4eSSatish Balay     }
1033ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
104e5c89e4eSSatish Balay   }
105e5c89e4eSSatish Balay   SWAP(vdx[0], vdx[right / 2], tmp);
106e5c89e4eSSatish Balay   vl   = v[vdx[0]];
107e5c89e4eSSatish Balay   last = 0;
108e5c89e4eSSatish Balay   for (i = 1; i <= right; i++) {
1099371c9d4SSatish Balay     if (v[vdx[i]] < vl) {
1109371c9d4SSatish Balay       last++;
1119371c9d4SSatish Balay       SWAP(vdx[last], vdx[i], tmp);
1129371c9d4SSatish Balay     }
113e5c89e4eSSatish Balay   }
114e5c89e4eSSatish Balay   SWAP(vdx[0], vdx[last], tmp);
1159566063dSJacob Faibussowitsch   PetscCall(PetscSortRealWithPermutation_Private(v, vdx, last - 1));
1169566063dSJacob Faibussowitsch   PetscCall(PetscSortRealWithPermutation_Private(v, vdx + last + 1, right - (last + 1)));
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118e5c89e4eSSatish Balay }
119e5c89e4eSSatish Balay 
120e5c89e4eSSatish Balay /*@
121811af0c4SBarry Smith   PetscSortRealWithPermutation - Computes the permutation of `PetscReal` that gives
122e5c89e4eSSatish Balay   a sorted sequence.
123e5c89e4eSSatish Balay 
124e5c89e4eSSatish Balay   Not Collective
125e5c89e4eSSatish Balay 
126e5c89e4eSSatish Balay   Input Parameters:
127e5c89e4eSSatish Balay + n   - number of values to sort
128e5c89e4eSSatish Balay . i   - values to sort
129*c435d973SPierre Jolivet - idx - permutation array. Must be initialized to 0:`n`-1 on input.
130e5c89e4eSSatish Balay 
131e5c89e4eSSatish Balay   Level: intermediate
132e5c89e4eSSatish Balay 
133811af0c4SBarry Smith   Note:
134*c435d973SPierre Jolivet   On output, `i` is unchanged and `idx[j]` is the position of the `j`th smallest `PetscReal` in `i`.
135e5c89e4eSSatish Balay 
136db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
137e5c89e4eSSatish Balay  @*/
PetscSortRealWithPermutation(PetscInt n,const PetscReal i[],PetscInt idx[])138d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRealWithPermutation(PetscInt n, const PetscReal i[], PetscInt idx[])
139d71ae5a4SJacob Faibussowitsch {
140e5c89e4eSSatish Balay   PetscInt  j, k, tmp;
141e5c89e4eSSatish Balay   PetscReal ik;
142e5c89e4eSSatish Balay 
143e5c89e4eSSatish Balay   PetscFunctionBegin;
144*c435d973SPierre Jolivet   if (n > 0) {
145*c435d973SPierre Jolivet     PetscAssertPointer(i, 2);
146*c435d973SPierre Jolivet     PetscAssertPointer(idx, 3);
147*c435d973SPierre Jolivet     PetscCheckIdentity(n, idx);
148*c435d973SPierre Jolivet   }
149e5c89e4eSSatish Balay   if (n < 8) {
150e5c89e4eSSatish Balay     for (k = 0; k < n; k++) {
151e5c89e4eSSatish Balay       ik = i[idx[k]];
152e5c89e4eSSatish Balay       for (j = k + 1; j < n; j++) {
153e5c89e4eSSatish Balay         if (ik > i[idx[j]]) {
154e5c89e4eSSatish Balay           SWAP(idx[k], idx[j], tmp);
155e5c89e4eSSatish Balay           ik = i[idx[k]];
156e5c89e4eSSatish Balay         }
157e5c89e4eSSatish Balay       }
158e5c89e4eSSatish Balay     }
159e5c89e4eSSatish Balay   } else {
1609566063dSJacob Faibussowitsch     PetscCall(PetscSortRealWithPermutation_Private(i, idx, n - 1));
161e5c89e4eSSatish Balay   }
1623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
163e5c89e4eSSatish Balay }
164e5c89e4eSSatish Balay 
PetscSortStrWithPermutation_Private(const char * v[],PetscInt vdx[],PetscInt right)165d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSortStrWithPermutation_Private(const char *v[], PetscInt vdx[], PetscInt right)
166d71ae5a4SJacob Faibussowitsch {
167e5c89e4eSSatish Balay   PetscInt    tmp, i, last;
168ace3abfcSBarry Smith   PetscBool   gt;
169e5c89e4eSSatish Balay   const char *vl;
170e5c89e4eSSatish Balay 
171e5c89e4eSSatish Balay   PetscFunctionBegin;
172e5c89e4eSSatish Balay   if (right <= 1) {
173e5c89e4eSSatish Balay     if (right == 1) {
1749566063dSJacob Faibussowitsch       PetscCall(PetscStrgrt(v[vdx[0]], v[vdx[1]], &gt));
175e5c89e4eSSatish Balay       if (gt) SWAP(vdx[0], vdx[1], tmp);
176e5c89e4eSSatish Balay     }
1773ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
178e5c89e4eSSatish Balay   }
179e5c89e4eSSatish Balay   SWAP(vdx[0], vdx[right / 2], tmp);
180e5c89e4eSSatish Balay   vl   = v[vdx[0]];
181e5c89e4eSSatish Balay   last = 0;
182e5c89e4eSSatish Balay   for (i = 1; i <= right; i++) {
1839566063dSJacob Faibussowitsch     PetscCall(PetscStrgrt(vl, v[vdx[i]], &gt));
1849371c9d4SSatish Balay     if (gt) {
1859371c9d4SSatish Balay       last++;
1869371c9d4SSatish Balay       SWAP(vdx[last], vdx[i], tmp);
1879371c9d4SSatish Balay     }
188e5c89e4eSSatish Balay   }
189e5c89e4eSSatish Balay   SWAP(vdx[0], vdx[last], tmp);
1909566063dSJacob Faibussowitsch   PetscCall(PetscSortStrWithPermutation_Private(v, vdx, last - 1));
1919566063dSJacob Faibussowitsch   PetscCall(PetscSortStrWithPermutation_Private(v, vdx + last + 1, right - (last + 1)));
1923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
193e5c89e4eSSatish Balay }
194e5c89e4eSSatish Balay 
195e5c89e4eSSatish Balay /*@C
196811af0c4SBarry Smith   PetscSortStrWithPermutation - Computes the permutation of strings that gives
197e5c89e4eSSatish Balay   a sorted sequence.
198e5c89e4eSSatish Balay 
199cc4c1da9SBarry Smith   Not Collective, No Fortran Support
200e5c89e4eSSatish Balay 
201e5c89e4eSSatish Balay   Input Parameters:
202e5c89e4eSSatish Balay + n   - number of values to sort
203e5c89e4eSSatish Balay . i   - values to sort
204*c435d973SPierre Jolivet - idx - permutation array. Must be initialized to 0:`n`-1 on input.
205e5c89e4eSSatish Balay 
206e5c89e4eSSatish Balay   Level: intermediate
207e5c89e4eSSatish Balay 
208811af0c4SBarry Smith   Note:
209*c435d973SPierre Jolivet   On output, `i` is unchanged and `idx[j]` is the position of the `j`th smallest `char *` in `i`.
210e5c89e4eSSatish Balay 
211db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()`
212e5c89e4eSSatish Balay  @*/
PetscSortStrWithPermutation(PetscInt n,const char * i[],PetscInt idx[])213d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortStrWithPermutation(PetscInt n, const char *i[], PetscInt idx[])
214d71ae5a4SJacob Faibussowitsch {
215e5c89e4eSSatish Balay   PetscInt    j, k, tmp;
216e5c89e4eSSatish Balay   const char *ik;
217ace3abfcSBarry Smith   PetscBool   gt;
218e5c89e4eSSatish Balay 
219e5c89e4eSSatish Balay   PetscFunctionBegin;
220*c435d973SPierre Jolivet   if (n > 0) {
221*c435d973SPierre Jolivet     PetscAssertPointer(i, 2);
222*c435d973SPierre Jolivet     PetscAssertPointer(idx, 3);
223*c435d973SPierre Jolivet     PetscCheckIdentity(n, idx);
224*c435d973SPierre Jolivet   }
225e5c89e4eSSatish Balay   if (n < 8) {
226e5c89e4eSSatish Balay     for (k = 0; k < n; k++) {
227e5c89e4eSSatish Balay       ik = i[idx[k]];
228e5c89e4eSSatish Balay       for (j = k + 1; j < n; j++) {
2299566063dSJacob Faibussowitsch         PetscCall(PetscStrgrt(ik, i[idx[j]], &gt));
230e5c89e4eSSatish Balay         if (gt) {
231e5c89e4eSSatish Balay           SWAP(idx[k], idx[j], tmp);
232e5c89e4eSSatish Balay           ik = i[idx[k]];
233e5c89e4eSSatish Balay         }
234e5c89e4eSSatish Balay       }
235e5c89e4eSSatish Balay     }
236e5c89e4eSSatish Balay   } else {
2379566063dSJacob Faibussowitsch     PetscCall(PetscSortStrWithPermutation_Private(i, idx, n - 1));
238e5c89e4eSSatish Balay   }
2393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
240e5c89e4eSSatish Balay }
241