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