Actual source code: sortip.c
1: /*
2: This file contains routines for sorting integers and doubles with a permutation array.
3: */
4: #include <petsc/private/petscimpl.h>
5: #include <petscsys.h>
7: #define SWAP(a, b, t) \
8: do { \
9: t = a; \
10: a = b; \
11: b = t; \
12: } while (0)
14: #if PetscDefined(USE_DEBUG)
15: #define PetscCheckIdentity(n, idx) \
16: do { \
17: 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: } while (0)
19: #else
20: #define PetscCheckIdentity(n, idx) (void)0
21: #endif
23: static PetscErrorCode PetscSortIntWithPermutation_Private(const PetscInt v[], PetscInt vdx[], PetscInt right)
24: {
25: PetscInt tmp, i, vl, last;
27: PetscFunctionBegin;
28: if (right <= 1) {
29: if (right == 1) {
30: if (v[vdx[0]] > v[vdx[1]]) SWAP(vdx[0], vdx[1], tmp);
31: }
32: PetscFunctionReturn(PETSC_SUCCESS);
33: }
34: SWAP(vdx[0], vdx[right / 2], tmp);
35: vl = v[vdx[0]];
36: last = 0;
37: for (i = 1; i <= right; i++) {
38: if (v[vdx[i]] < vl) {
39: last++;
40: SWAP(vdx[last], vdx[i], tmp);
41: }
42: }
43: SWAP(vdx[0], vdx[last], tmp);
44: PetscCall(PetscSortIntWithPermutation_Private(v, vdx, last - 1));
45: PetscCall(PetscSortIntWithPermutation_Private(v, vdx + last + 1, right - (last + 1)));
46: PetscFunctionReturn(PETSC_SUCCESS);
47: }
49: /*@
50: PetscSortIntWithPermutation - Computes the permutation of `PetscInt` that gives
51: a sorted sequence.
53: Not Collective
55: Input Parameters:
56: + n - number of values to sort
57: . i - values to sort
58: - idx - permutation array. Must be initialized to 0:`n`-1 on input.
60: Level: intermediate
62: Note:
63: On output, `i` is unchanged and `idx[j]` is the position of the `j`th smallest `PetscInt` in `i`.
65: .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()`, `PetscSortIntWithArray()`
66: @*/
67: PetscErrorCode PetscSortIntWithPermutation(PetscInt n, const PetscInt i[], PetscInt idx[])
68: {
69: PetscInt j, k, tmp, ik;
71: PetscFunctionBegin;
72: if (n > 0) {
73: PetscAssertPointer(i, 2);
74: PetscAssertPointer(idx, 3);
75: PetscCheckIdentity(n, idx);
76: }
77: if (n < 8) {
78: for (k = 0; k < n; k++) {
79: ik = i[idx[k]];
80: for (j = k + 1; j < n; j++) {
81: if (ik > i[idx[j]]) {
82: SWAP(idx[k], idx[j], tmp);
83: ik = i[idx[k]];
84: }
85: }
86: }
87: } else {
88: PetscCall(PetscSortIntWithPermutation_Private(i, idx, n - 1));
89: }
90: PetscFunctionReturn(PETSC_SUCCESS);
91: }
93: static PetscErrorCode PetscSortRealWithPermutation_Private(const PetscReal v[], PetscInt vdx[], PetscInt right)
94: {
95: PetscReal vl;
96: PetscInt tmp, i, last;
98: PetscFunctionBegin;
99: if (right <= 1) {
100: if (right == 1) {
101: if (v[vdx[0]] > v[vdx[1]]) SWAP(vdx[0], vdx[1], tmp);
102: }
103: PetscFunctionReturn(PETSC_SUCCESS);
104: }
105: SWAP(vdx[0], vdx[right / 2], tmp);
106: vl = v[vdx[0]];
107: last = 0;
108: for (i = 1; i <= right; i++) {
109: if (v[vdx[i]] < vl) {
110: last++;
111: SWAP(vdx[last], vdx[i], tmp);
112: }
113: }
114: SWAP(vdx[0], vdx[last], tmp);
115: PetscCall(PetscSortRealWithPermutation_Private(v, vdx, last - 1));
116: PetscCall(PetscSortRealWithPermutation_Private(v, vdx + last + 1, right - (last + 1)));
117: PetscFunctionReturn(PETSC_SUCCESS);
118: }
120: /*@
121: PetscSortRealWithPermutation - Computes the permutation of `PetscReal` that gives
122: a sorted sequence.
124: Not Collective
126: Input Parameters:
127: + n - number of values to sort
128: . i - values to sort
129: - idx - permutation array. Must be initialized to 0:`n`-1 on input.
131: Level: intermediate
133: Note:
134: On output, `i` is unchanged and `idx[j]` is the position of the `j`th smallest `PetscReal` in `i`.
136: .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
137: @*/
138: PetscErrorCode PetscSortRealWithPermutation(PetscInt n, const PetscReal i[], PetscInt idx[])
139: {
140: PetscInt j, k, tmp;
141: PetscReal ik;
143: PetscFunctionBegin;
144: if (n > 0) {
145: PetscAssertPointer(i, 2);
146: PetscAssertPointer(idx, 3);
147: PetscCheckIdentity(n, idx);
148: }
149: if (n < 8) {
150: for (k = 0; k < n; k++) {
151: ik = i[idx[k]];
152: for (j = k + 1; j < n; j++) {
153: if (ik > i[idx[j]]) {
154: SWAP(idx[k], idx[j], tmp);
155: ik = i[idx[k]];
156: }
157: }
158: }
159: } else {
160: PetscCall(PetscSortRealWithPermutation_Private(i, idx, n - 1));
161: }
162: PetscFunctionReturn(PETSC_SUCCESS);
163: }
165: static PetscErrorCode PetscSortStrWithPermutation_Private(const char *v[], PetscInt vdx[], PetscInt right)
166: {
167: PetscInt tmp, i, last;
168: PetscBool gt;
169: const char *vl;
171: PetscFunctionBegin;
172: if (right <= 1) {
173: if (right == 1) {
174: PetscCall(PetscStrgrt(v[vdx[0]], v[vdx[1]], >));
175: if (gt) SWAP(vdx[0], vdx[1], tmp);
176: }
177: PetscFunctionReturn(PETSC_SUCCESS);
178: }
179: SWAP(vdx[0], vdx[right / 2], tmp);
180: vl = v[vdx[0]];
181: last = 0;
182: for (i = 1; i <= right; i++) {
183: PetscCall(PetscStrgrt(vl, v[vdx[i]], >));
184: if (gt) {
185: last++;
186: SWAP(vdx[last], vdx[i], tmp);
187: }
188: }
189: SWAP(vdx[0], vdx[last], tmp);
190: PetscCall(PetscSortStrWithPermutation_Private(v, vdx, last - 1));
191: PetscCall(PetscSortStrWithPermutation_Private(v, vdx + last + 1, right - (last + 1)));
192: PetscFunctionReturn(PETSC_SUCCESS);
193: }
195: /*@C
196: PetscSortStrWithPermutation - Computes the permutation of strings that gives
197: a sorted sequence.
199: Not Collective
201: Input Parameters:
202: + n - number of values to sort
203: . i - values to sort
204: - idx - permutation array. Must be initialized to 0:`n`-1 on input.
206: Level: intermediate
208: Note:
209: On output, `i` is unchanged and `idx[j]` is the position of the `j`th smallest `char *` in `i`.
211: .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()`
212: @*/
213: PetscErrorCode PetscSortStrWithPermutation(PetscInt n, const char *i[], PetscInt idx[])
214: {
215: PetscInt j, k, tmp;
216: const char *ik;
217: PetscBool gt;
219: PetscFunctionBegin;
220: if (n > 0) {
221: PetscAssertPointer(i, 2);
222: PetscAssertPointer(idx, 3);
223: PetscCheckIdentity(n, idx);
224: }
225: if (n < 8) {
226: for (k = 0; k < n; k++) {
227: ik = i[idx[k]];
228: for (j = k + 1; j < n; j++) {
229: PetscCall(PetscStrgrt(ik, i[idx[j]], >));
230: if (gt) {
231: SWAP(idx[k], idx[j], tmp);
232: ik = i[idx[k]];
233: }
234: }
235: }
236: } else {
237: PetscCall(PetscSortStrWithPermutation_Private(i, idx, n - 1));
238: }
239: PetscFunctionReturn(PETSC_SUCCESS);
240: }