Actual source code: sortip.c

  1: /*
  2:    This file contains routines for sorting integers and doubles with a permutation array.

  4:    The word "register"  in this code is used to identify data that is not
  5:    aliased.  For some compilers, this can cause the compiler to fail to
  6:    place inner-loop variables into registers.
  7:  */
  8: #include <petscsys.h>

 10: #define SWAP(a, b, t) \
 11:   do { \
 12:     t = a; \
 13:     a = b; \
 14:     b = t; \
 15:   } while (0)

 17: static PetscErrorCode PetscSortIntWithPermutation_Private(const PetscInt v[], PetscInt vdx[], PetscInt right)
 18: {
 19:   PetscInt tmp, i, vl, last;

 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(PETSC_SUCCESS);
 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(PETSC_SUCCESS);
 41: }

 43: /*@
 44:   PetscSortIntWithPermutation - Computes the permutation of `PetscInt` that gives
 45:   a sorted sequence.

 47:   Not Collective

 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.

 54:   Level: intermediate

 56:   Note:
 57:   On output i is unchanged and idx[i] is the position of the i-th smallest index in i.

 59: .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()`, `PetscSortIntWithArray()`
 60:  @*/
 61: PetscErrorCode PetscSortIntWithPermutation(PetscInt n, const PetscInt i[], PetscInt idx[])
 62: {
 63:   PetscInt j, k, tmp, ik;

 65:   PetscFunctionBegin;
 66:   if (n < 8) {
 67:     for (k = 0; k < n; k++) {
 68:       ik = i[idx[k]];
 69:       for (j = k + 1; j < n; j++) {
 70:         if (ik > i[idx[j]]) {
 71:           SWAP(idx[k], idx[j], tmp);
 72:           ik = i[idx[k]];
 73:         }
 74:       }
 75:     }
 76:   } else {
 77:     PetscCall(PetscSortIntWithPermutation_Private(i, idx, n - 1));
 78:   }
 79:   PetscFunctionReturn(PETSC_SUCCESS);
 80: }

 82: /* ---------------------------------------------------------------------- */

 84: static PetscErrorCode PetscSortRealWithPermutation_Private(const PetscReal v[], PetscInt vdx[], PetscInt right)
 85: {
 86:   PetscReal vl;
 87:   PetscInt  tmp, i, last;

 89:   PetscFunctionBegin;
 90:   if (right <= 1) {
 91:     if (right == 1) {
 92:       if (v[vdx[0]] > v[vdx[1]]) SWAP(vdx[0], vdx[1], tmp);
 93:     }
 94:     PetscFunctionReturn(PETSC_SUCCESS);
 95:   }
 96:   SWAP(vdx[0], vdx[right / 2], tmp);
 97:   vl   = v[vdx[0]];
 98:   last = 0;
 99:   for (i = 1; i <= right; i++) {
100:     if (v[vdx[i]] < vl) {
101:       last++;
102:       SWAP(vdx[last], vdx[i], tmp);
103:     }
104:   }
105:   SWAP(vdx[0], vdx[last], tmp);
106:   PetscCall(PetscSortRealWithPermutation_Private(v, vdx, last - 1));
107:   PetscCall(PetscSortRealWithPermutation_Private(v, vdx + last + 1, right - (last + 1)));
108:   PetscFunctionReturn(PETSC_SUCCESS);
109: }

111: /*@
112:   PetscSortRealWithPermutation - Computes the permutation of `PetscReal` that gives
113:   a sorted sequence.

115:   Not Collective

117:   Input Parameters:
118: + n   - number of values to sort
119: . i   - values to sort
120: - idx - permutation array.  Must be initialized to 0:n-1 on input.

122:   Level: intermediate

124:   Note:
125:   i is unchanged on output.

127: .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
128:  @*/
129: PetscErrorCode PetscSortRealWithPermutation(PetscInt n, const PetscReal i[], PetscInt idx[])
130: {
131:   PetscInt  j, k, tmp;
132:   PetscReal ik;

134:   PetscFunctionBegin;
135:   if (n < 8) {
136:     for (k = 0; k < n; k++) {
137:       ik = i[idx[k]];
138:       for (j = k + 1; j < n; j++) {
139:         if (ik > i[idx[j]]) {
140:           SWAP(idx[k], idx[j], tmp);
141:           ik = i[idx[k]];
142:         }
143:       }
144:     }
145:   } else {
146:     PetscCall(PetscSortRealWithPermutation_Private(i, idx, n - 1));
147:   }
148:   PetscFunctionReturn(PETSC_SUCCESS);
149: }

151: static PetscErrorCode PetscSortStrWithPermutation_Private(const char *v[], PetscInt vdx[], PetscInt right)
152: {
153:   PetscInt    tmp, i, last;
154:   PetscBool   gt;
155:   const char *vl;

157:   PetscFunctionBegin;
158:   if (right <= 1) {
159:     if (right == 1) {
160:       PetscCall(PetscStrgrt(v[vdx[0]], v[vdx[1]], &gt));
161:       if (gt) SWAP(vdx[0], vdx[1], tmp);
162:     }
163:     PetscFunctionReturn(PETSC_SUCCESS);
164:   }
165:   SWAP(vdx[0], vdx[right / 2], tmp);
166:   vl   = v[vdx[0]];
167:   last = 0;
168:   for (i = 1; i <= right; i++) {
169:     PetscCall(PetscStrgrt(vl, v[vdx[i]], &gt));
170:     if (gt) {
171:       last++;
172:       SWAP(vdx[last], vdx[i], tmp);
173:     }
174:   }
175:   SWAP(vdx[0], vdx[last], tmp);
176:   PetscCall(PetscSortStrWithPermutation_Private(v, vdx, last - 1));
177:   PetscCall(PetscSortStrWithPermutation_Private(v, vdx + last + 1, right - (last + 1)));
178:   PetscFunctionReturn(PETSC_SUCCESS);
179: }

181: /*@C
182:   PetscSortStrWithPermutation - Computes the permutation of strings that gives
183:   a sorted sequence.

185:   Not Collective

187:   Input Parameters:
188: + n   - number of values to sort
189: . i   - values to sort
190: - idx - permutation array.  Must be initialized to 0:n-1 on input.

192:   Level: intermediate

194:   Note:
195:   i is unchanged on output.

197: .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()`
198:  @*/
199: PetscErrorCode PetscSortStrWithPermutation(PetscInt n, const char *i[], PetscInt idx[])
200: {
201:   PetscInt    j, k, tmp;
202:   const char *ik;
203:   PetscBool   gt;

205:   PetscFunctionBegin;
206:   if (n < 8) {
207:     for (k = 0; k < n; k++) {
208:       ik = i[idx[k]];
209:       for (j = k + 1; j < n; j++) {
210:         PetscCall(PetscStrgrt(ik, i[idx[j]], &gt));
211:         if (gt) {
212:           SWAP(idx[k], idx[j], tmp);
213:           ik = i[idx[k]];
214:         }
215:       }
216:     }
217:   } else {
218:     PetscCall(PetscSortStrWithPermutation_Private(i, idx, n - 1));
219:   }
220:   PetscFunctionReturn(PETSC_SUCCESS);
221: }