xref: /petsc/src/sys/utils/sortip.c (revision 21e3ffae2f3b73c0bd738cf6d0a809700fc04bb0)
1 
2 /*
3    This file contains routines for sorting integers and doubles with a permutation array.
4 
5    The word "register"  in this code is used to identify data that is not
6    aliased.  For some compilers, this can cause the compiler to fail to
7    place inner-loop variables into registers.
8  */
9 #include <petscsys.h> /*I  "petscsys.h"  I*/
10 
11 #define SWAP(a, b, t) \
12   { \
13     t = a; \
14     a = b; \
15     b = t; \
16   }
17 
18 static PetscErrorCode PetscSortIntWithPermutation_Private(const PetscInt v[], PetscInt vdx[], PetscInt right)
19 {
20   PetscInt tmp, i, vl, last;
21 
22   PetscFunctionBegin;
23   if (right <= 1) {
24     if (right == 1) {
25       if (v[vdx[0]] > v[vdx[1]]) SWAP(vdx[0], vdx[1], tmp);
26     }
27     PetscFunctionReturn(PETSC_SUCCESS);
28   }
29   SWAP(vdx[0], vdx[right / 2], tmp);
30   vl   = v[vdx[0]];
31   last = 0;
32   for (i = 1; i <= right; i++) {
33     if (v[vdx[i]] < vl) {
34       last++;
35       SWAP(vdx[last], vdx[i], tmp);
36     }
37   }
38   SWAP(vdx[0], vdx[last], tmp);
39   PetscCall(PetscSortIntWithPermutation_Private(v, vdx, last - 1));
40   PetscCall(PetscSortIntWithPermutation_Private(v, vdx + last + 1, right - (last + 1)));
41   PetscFunctionReturn(PETSC_SUCCESS);
42 }
43 
44 /*@
45    PetscSortIntWithPermutation - Computes the permutation of `PetscInt` that gives
46    a sorted sequence.
47 
48    Not Collective
49 
50    Input Parameters:
51 +  n  - number of values to sort
52 .  i  - values to sort
53 -  idx - permutation array.  Must be initialized to 0:n-1 on input.
54 
55    Level: intermediate
56 
57    Note:
58    On output i is unchanged and idx[i] is the position of the i-th smallest index in i.
59 
60 .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()`, `PetscSortIntWithArray()`
61  @*/
62 PetscErrorCode PetscSortIntWithPermutation(PetscInt n, const PetscInt i[], PetscInt idx[])
63 {
64   PetscInt j, k, tmp, ik;
65 
66   PetscFunctionBegin;
67   if (n < 8) {
68     for (k = 0; k < n; k++) {
69       ik = i[idx[k]];
70       for (j = k + 1; j < n; j++) {
71         if (ik > i[idx[j]]) {
72           SWAP(idx[k], idx[j], tmp);
73           ik = i[idx[k]];
74         }
75       }
76     }
77   } else {
78     PetscCall(PetscSortIntWithPermutation_Private(i, idx, n - 1));
79   }
80   PetscFunctionReturn(PETSC_SUCCESS);
81 }
82 
83 /* ---------------------------------------------------------------------- */
84 
85 static PetscErrorCode PetscSortRealWithPermutation_Private(const PetscReal v[], PetscInt vdx[], PetscInt right)
86 {
87   PetscReal vl;
88   PetscInt  tmp, i, last;
89 
90   PetscFunctionBegin;
91   if (right <= 1) {
92     if (right == 1) {
93       if (v[vdx[0]] > v[vdx[1]]) SWAP(vdx[0], vdx[1], tmp);
94     }
95     PetscFunctionReturn(PETSC_SUCCESS);
96   }
97   SWAP(vdx[0], vdx[right / 2], tmp);
98   vl   = v[vdx[0]];
99   last = 0;
100   for (i = 1; i <= right; i++) {
101     if (v[vdx[i]] < vl) {
102       last++;
103       SWAP(vdx[last], vdx[i], tmp);
104     }
105   }
106   SWAP(vdx[0], vdx[last], tmp);
107   PetscCall(PetscSortRealWithPermutation_Private(v, vdx, last - 1));
108   PetscCall(PetscSortRealWithPermutation_Private(v, vdx + last + 1, right - (last + 1)));
109   PetscFunctionReturn(PETSC_SUCCESS);
110 }
111 
112 /*@
113    PetscSortRealWithPermutation - Computes the permutation of `PetscReal` that gives
114    a sorted sequence.
115 
116    Not Collective
117 
118    Input Parameters:
119 +  n  - number of values to sort
120 .  i  - values to sort
121 -  idx - permutation array.  Must be initialized to 0:n-1 on input.
122 
123    Level: intermediate
124 
125    Note:
126    i is unchanged on output.
127 
128 .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
129  @*/
130 PetscErrorCode PetscSortRealWithPermutation(PetscInt n, const PetscReal i[], PetscInt idx[])
131 {
132   PetscInt  j, k, tmp;
133   PetscReal ik;
134 
135   PetscFunctionBegin;
136   if (n < 8) {
137     for (k = 0; k < n; k++) {
138       ik = i[idx[k]];
139       for (j = k + 1; j < n; j++) {
140         if (ik > i[idx[j]]) {
141           SWAP(idx[k], idx[j], tmp);
142           ik = i[idx[k]];
143         }
144       }
145     }
146   } else {
147     PetscCall(PetscSortRealWithPermutation_Private(i, idx, n - 1));
148   }
149   PetscFunctionReturn(PETSC_SUCCESS);
150 }
151 
152 static PetscErrorCode PetscSortStrWithPermutation_Private(const char *v[], PetscInt vdx[], PetscInt right)
153 {
154   PetscInt    tmp, i, last;
155   PetscBool   gt;
156   const char *vl;
157 
158   PetscFunctionBegin;
159   if (right <= 1) {
160     if (right == 1) {
161       PetscCall(PetscStrgrt(v[vdx[0]], v[vdx[1]], &gt));
162       if (gt) SWAP(vdx[0], vdx[1], tmp);
163     }
164     PetscFunctionReturn(PETSC_SUCCESS);
165   }
166   SWAP(vdx[0], vdx[right / 2], tmp);
167   vl   = v[vdx[0]];
168   last = 0;
169   for (i = 1; i <= right; i++) {
170     PetscCall(PetscStrgrt(vl, v[vdx[i]], &gt));
171     if (gt) {
172       last++;
173       SWAP(vdx[last], vdx[i], tmp);
174     }
175   }
176   SWAP(vdx[0], vdx[last], tmp);
177   PetscCall(PetscSortStrWithPermutation_Private(v, vdx, last - 1));
178   PetscCall(PetscSortStrWithPermutation_Private(v, vdx + last + 1, right - (last + 1)));
179   PetscFunctionReturn(PETSC_SUCCESS);
180 }
181 
182 /*@C
183    PetscSortStrWithPermutation - Computes the permutation of strings that gives
184    a sorted sequence.
185 
186    Not Collective
187 
188    Input Parameters:
189 +  n  - number of values to sort
190 .  i  - values to sort
191 -  idx - permutation array.  Must be initialized to 0:n-1 on input.
192 
193    Level: intermediate
194 
195    Note:
196    i is unchanged on output.
197 
198 .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()`
199  @*/
200 PetscErrorCode PetscSortStrWithPermutation(PetscInt n, const char *i[], PetscInt idx[])
201 {
202   PetscInt    j, k, tmp;
203   const char *ik;
204   PetscBool   gt;
205 
206   PetscFunctionBegin;
207   if (n < 8) {
208     for (k = 0; k < n; k++) {
209       ik = i[idx[k]];
210       for (j = k + 1; j < n; j++) {
211         PetscCall(PetscStrgrt(ik, i[idx[j]], &gt));
212         if (gt) {
213           SWAP(idx[k], idx[j], tmp);
214           ik = i[idx[k]];
215         }
216       }
217     }
218   } else {
219     PetscCall(PetscSortStrWithPermutation_Private(i, idx, n - 1));
220   }
221   PetscFunctionReturn(PETSC_SUCCESS);
222 }
223