Actual source code: sorti.c

  1: /*
  2:    This file contains routines for sorting integers. Values are sorted in place.
  3:    One can use src/sys/tests/ex52.c for benchmarking.
  4:  */
  5: #include <petsc/private/petscimpl.h>
  6: #include <petsc/private/hashseti.h>

  8: #define MEDIAN3(v, a, b, c) (v[a] < v[b] ? (v[b] < v[c] ? (b) : (v[a] < v[c] ? (c) : (a))) : (v[c] < v[b] ? (b) : (v[a] < v[c] ? (a) : (c))))

 10: #define MEDIAN(v, right) MEDIAN3(v, right / 4, right / 2, right / 4 * 3)

 12: /* Swap one, two or three pairs. Each pair can have its own type */
 13: #define SWAP1(a, b, t1) \
 14:   do { \
 15:     t1 = a; \
 16:     a  = b; \
 17:     b  = t1; \
 18:   } while (0)
 19: #define SWAP2(a, b, c, d, t1, t2) \
 20:   do { \
 21:     t1 = a; \
 22:     a  = b; \
 23:     b  = t1; \
 24:     t2 = c; \
 25:     c  = d; \
 26:     d  = t2; \
 27:   } while (0)
 28: #define SWAP3(a, b, c, d, e, f, t1, t2, t3) \
 29:   do { \
 30:     t1 = a; \
 31:     a  = b; \
 32:     b  = t1; \
 33:     t2 = c; \
 34:     c  = d; \
 35:     d  = t2; \
 36:     t3 = e; \
 37:     e  = f; \
 38:     f  = t3; \
 39:   } while (0)

 41: /* Swap a & b, *c & *d. c, d, t2 are pointers to a type of size <siz> */
 42: #define SWAP2Data(a, b, c, d, t1, t2, siz) \
 43:   do { \
 44:     t1 = a; \
 45:     a  = b; \
 46:     b  = t1; \
 47:     PetscCall(PetscMemcpy(t2, c, siz)); \
 48:     PetscCall(PetscMemcpy(c, d, siz)); \
 49:     PetscCall(PetscMemcpy(d, t2, siz)); \
 50:   } while (0)

 52: /*
 53:    Partition X[lo,hi] into two parts: X[lo,l) <= pivot; X[r,hi] > pivot

 55:    Input Parameters:
 56:     + X         - array to partition
 57:     . pivot     - a pivot of X[]
 58:     . t1        - temp variable for X
 59:     - lo,hi     - lower and upper bound of the array

 61:    Output Parameters:
 62:     + l,r       - of type PetscInt

 64:    Note:
 65:     The TwoWayPartition2/3 variants also partition other arrays along with X.
 66:     These arrays can have different types, so they provide their own temp t2,t3
 67:  */
 68: #define TwoWayPartition1(X, pivot, t1, lo, hi, l, r) \
 69:   do { \
 70:     l = lo; \
 71:     r = hi; \
 72:     while (1) { \
 73:       while (X[l] < pivot) l++; \
 74:       while (X[r] > pivot) r--; \
 75:       if (l >= r) { \
 76:         r++; \
 77:         break; \
 78:       } \
 79:       SWAP1(X[l], X[r], t1); \
 80:       l++; \
 81:       r--; \
 82:     } \
 83:   } while (0)

 85: /*
 86:    Partition X[lo,hi] into two parts: X[lo,l) >= pivot; X[r,hi] < pivot

 88:    Input Parameters:
 89:     + X         - array to partition
 90:     . pivot     - a pivot of X[]
 91:     . t1        - temp variable for X
 92:     - lo,hi     - lower and upper bound of the array

 94:    Output Parameters:
 95:     + l,r       - of type PetscInt

 97:    Note:
 98:     The TwoWayPartition2/3 variants also partition other arrays along with X.
 99:     These arrays can have different types, so they provide their own temp t2,t3
100:  */
101: #define TwoWayPartitionReverse1(X, pivot, t1, lo, hi, l, r) \
102:   do { \
103:     l = lo; \
104:     r = hi; \
105:     while (1) { \
106:       while (X[l] > pivot) l++; \
107:       while (X[r] < pivot) r--; \
108:       if (l >= r) { \
109:         r++; \
110:         break; \
111:       } \
112:       SWAP1(X[l], X[r], t1); \
113:       l++; \
114:       r--; \
115:     } \
116:   } while (0)

118: #define TwoWayPartition2(X, Y, pivot, t1, t2, lo, hi, l, r) \
119:   do { \
120:     l = lo; \
121:     r = hi; \
122:     while (1) { \
123:       while (X[l] < pivot) l++; \
124:       while (X[r] > pivot) r--; \
125:       if (l >= r) { \
126:         r++; \
127:         break; \
128:       } \
129:       SWAP2(X[l], X[r], Y[l], Y[r], t1, t2); \
130:       l++; \
131:       r--; \
132:     } \
133:   } while (0)

135: #define TwoWayPartition3(X, Y, Z, pivot, t1, t2, t3, lo, hi, l, r) \
136:   do { \
137:     l = lo; \
138:     r = hi; \
139:     while (1) { \
140:       while (X[l] < pivot) l++; \
141:       while (X[r] > pivot) r--; \
142:       if (l >= r) { \
143:         r++; \
144:         break; \
145:       } \
146:       SWAP3(X[l], X[r], Y[l], Y[r], Z[l], Z[r], t1, t2, t3); \
147:       l++; \
148:       r--; \
149:     } \
150:   } while (0)

152: /* Templates for similar functions used below */
153: #define QuickSort1(FuncName, X, n, pivot, t1) \
154:   do { \
155:     PetscCount i, j, p, l, r, hi = n - 1; \
156:     if (n < 8) { \
157:       for (i = 0; i < n; i++) { \
158:         pivot = X[i]; \
159:         for (j = i + 1; j < n; j++) { \
160:           if (pivot > X[j]) { \
161:             SWAP1(X[i], X[j], t1); \
162:             pivot = X[i]; \
163:           } \
164:         } \
165:       } \
166:     } else { \
167:       p     = MEDIAN(X, hi); \
168:       pivot = X[p]; \
169:       TwoWayPartition1(X, pivot, t1, 0, hi, l, r); \
170:       PetscCall(FuncName(l, X)); \
171:       PetscCall(FuncName(hi - r + 1, X + r)); \
172:     } \
173:   } while (0)

175: /* Templates for similar functions used below */
176: #define QuickSortReverse1(FuncName, X, n, pivot, t1) \
177:   do { \
178:     PetscCount i, j, p, l, r, hi = n - 1; \
179:     if (n < 8) { \
180:       for (i = 0; i < n; i++) { \
181:         pivot = X[i]; \
182:         for (j = i + 1; j < n; j++) { \
183:           if (pivot < X[j]) { \
184:             SWAP1(X[i], X[j], t1); \
185:             pivot = X[i]; \
186:           } \
187:         } \
188:       } \
189:     } else { \
190:       p     = MEDIAN(X, hi); \
191:       pivot = X[p]; \
192:       TwoWayPartitionReverse1(X, pivot, t1, 0, hi, l, r); \
193:       PetscCall(FuncName(l, X)); \
194:       PetscCall(FuncName(hi - r + 1, X + r)); \
195:     } \
196:   } while (0)

198: #define QuickSort2(FuncName, X, Y, n, pivot, t1, t2) \
199:   do { \
200:     PetscCount i, j, p, l, r, hi = n - 1; \
201:     if (n < 8) { \
202:       for (i = 0; i < n; i++) { \
203:         pivot = X[i]; \
204:         for (j = i + 1; j < n; j++) { \
205:           if (pivot > X[j]) { \
206:             SWAP2(X[i], X[j], Y[i], Y[j], t1, t2); \
207:             pivot = X[i]; \
208:           } \
209:         } \
210:       } \
211:     } else { \
212:       p     = MEDIAN(X, hi); \
213:       pivot = X[p]; \
214:       TwoWayPartition2(X, Y, pivot, t1, t2, 0, hi, l, r); \
215:       PetscCall(FuncName(l, X, Y)); \
216:       PetscCall(FuncName(hi - r + 1, X + r, Y + r)); \
217:     } \
218:   } while (0)

220: #define QuickSort3(FuncName, X, Y, Z, n, pivot, t1, t2, t3) \
221:   do { \
222:     PetscCount i, j, p, l, r, hi = n - 1; \
223:     if (n < 8) { \
224:       for (i = 0; i < n; i++) { \
225:         pivot = X[i]; \
226:         for (j = i + 1; j < n; j++) { \
227:           if (pivot > X[j]) { \
228:             SWAP3(X[i], X[j], Y[i], Y[j], Z[i], Z[j], t1, t2, t3); \
229:             pivot = X[i]; \
230:           } \
231:         } \
232:       } \
233:     } else { \
234:       p     = MEDIAN(X, hi); \
235:       pivot = X[p]; \
236:       TwoWayPartition3(X, Y, Z, pivot, t1, t2, t3, 0, hi, l, r); \
237:       PetscCall(FuncName(l, X, Y, Z)); \
238:       PetscCall(FuncName(hi - r + 1, X + r, Y + r, Z + r)); \
239:     } \
240:   } while (0)

242: /*@
243:   PetscSortedInt - Determines whether the `PetscInt` array is sorted.

245:   Not Collective

247:   Input Parameters:
248: + n - number of values
249: - X - array of `PetscInt`

251:   Output Parameter:
252: . sorted - flag whether the array is sorted

254:   Level: intermediate

256: .seealso: `PetscSortInt()`, `PetscSortedMPIInt()`, `PetscSortedReal()`
257: @*/
258: PetscErrorCode PetscSortedInt(PetscCount n, const PetscInt X[], PetscBool *sorted)
259: {
260:   PetscFunctionBegin;
261:   if (n) PetscAssertPointer(X, 2);
262:   PetscAssertPointer(sorted, 3);
263:   PetscSorted(n, X, *sorted);
264:   PetscFunctionReturn(PETSC_SUCCESS);
265: }

267: /*@
268:   PetscSortedInt64 - Determines whether the `PetscInt64` array is sorted.

270:   Not Collective

272:   Input Parameters:
273: + n - number of values
274: - X - array of `PetscInt64`

276:   Output Parameter:
277: . sorted - flag whether the array is sorted

279:   Level: intermediate

281: .seealso: `PetscSortInt64()`, `PetscSortInt()`, `PetscSortedMPIInt()`, `PetscSortedReal()`
282: @*/
283: PetscErrorCode PetscSortedInt64(PetscCount n, const PetscInt64 X[], PetscBool *sorted)
284: {
285:   PetscFunctionBegin;
286:   if (n) PetscAssertPointer(X, 2);
287:   PetscAssertPointer(sorted, 3);
288:   PetscSorted(n, X, *sorted);
289:   PetscFunctionReturn(PETSC_SUCCESS);
290: }

292: /*@
293:   PetscSortInt - Sorts an array of `PetscInt` in place in increasing order.

295:   Not Collective

297:   Input Parameters:
298: + n - number of values
299: - X - array of `PetscInt`

301:   Note:
302:   This function serves as an alternative to `PetscIntSortSemiOrdered()`, and may perform faster especially if the array
303:   is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
304:   code to see which routine is fastest.

306:   Level: intermediate

308: .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
309: @*/
310: PetscErrorCode PetscSortInt(PetscCount n, PetscInt X[])
311: {
312:   PetscInt pivot, t1, N;

314:   PetscFunctionBegin;
315:   if (n) PetscAssertPointer(X, 2);
316:   PetscCall(PetscIntCast(n, &N));
317:   QuickSort1(PetscSortInt, X, N, pivot, t1);
318:   PetscFunctionReturn(PETSC_SUCCESS);
319: }

321: /*@
322:   PetscSortInt64 - Sorts an array of `PetscInt64` in place in increasing order.

324:   Not Collective

326:   Input Parameters:
327: + n - number of values
328: - X - array of `PetscInt64`

330:   Notes:
331:   This function sorts `PetscInt64`s assumed to be in completely random order

333:   Level: intermediate

335: .seealso: `PetscSortInt()`
336: @*/
337: PetscErrorCode PetscSortInt64(PetscCount n, PetscInt64 X[])
338: {
339:   PetscCount pivot, t1;

341:   PetscFunctionBegin;
342:   if (n) PetscAssertPointer(X, 2);
343:   QuickSort1(PetscSortInt64, X, n, pivot, t1);
344:   PetscFunctionReturn(PETSC_SUCCESS);
345: }

347: /*@
348:   PetscSortCount - Sorts an array of `PetscCount` in place in increasing order.

350:   Not Collective

352:   Input Parameters:
353: + n - number of values
354: - X - array of `PetscCount`

356:   Notes:
357:   This function sorts `PetscCount`s assumed to be in completely random order

359:   Level: intermediate

361: .seealso: `PetscSortInt()`
362: @*/
363: PetscErrorCode PetscSortCount(PetscCount n, PetscCount X[])
364: {
365:   PetscCount pivot, t1;

367:   PetscFunctionBegin;
368:   if (n) PetscAssertPointer(X, 2);
369:   QuickSort1(PetscSortCount, X, n, pivot, t1);
370:   PetscFunctionReturn(PETSC_SUCCESS);
371: }

373: /*@
374:   PetscSortReverseInt - Sorts an array of `PetscInt` in place in decreasing order.

376:   Not Collective

378:   Input Parameters:
379: + n - number of values
380: - X - array of `PetscInt`

382:   Level: intermediate

384: .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithPermutation()`
385: @*/
386: PetscErrorCode PetscSortReverseInt(PetscCount n, PetscInt X[])
387: {
388:   PetscInt pivot, t1;

390:   PetscFunctionBegin;
391:   if (n) PetscAssertPointer(X, 2);
392:   QuickSortReverse1(PetscSortReverseInt, X, n, pivot, t1);
393:   PetscFunctionReturn(PETSC_SUCCESS);
394: }

396: /*@
397:   PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted `PetscInt` array

399:   Not Collective

401:   Input Parameters:
402: + n - number of values
403: - X - sorted array of `PetscInt`

405:   Output Parameter:
406: . n - number of non-redundant values

408:   Level: intermediate

410: .seealso: `PetscSortInt()`
411: @*/
412: PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n, PetscInt X[])
413: {
414:   PetscInt i, s = 0, N = *n, b = 0;

416:   PetscFunctionBegin;
417:   PetscAssertPointer(n, 1);
418:   PetscCheckSorted(*n, X);
419:   for (i = 0; i < N - 1; i++) {
420:     if (X[b + s + 1] != X[b]) {
421:       X[b + 1] = X[b + s + 1];
422:       b++;
423:     } else s++;
424:   }
425:   *n = N - s;
426:   PetscFunctionReturn(PETSC_SUCCESS);
427: }

429: /*@
430:   PetscSortedCheckDupsInt - Checks if a sorted `PetscInt` array has duplicates

432:   Not Collective

434:   Input Parameters:
435: + n - number of values
436: - X - sorted array of `PetscInt`

438:   Output Parameter:
439: . flg - True if the array has duplications, otherwise false

441:   Level: intermediate

443: .seealso: `PetscSortInt()`, `PetscCheckDupsInt()`, `PetscSortRemoveDupsInt()`, `PetscSortedRemoveDupsInt()`
444: @*/
445: PetscErrorCode PetscSortedCheckDupsInt(PetscCount n, const PetscInt X[], PetscBool *flg)
446: {
447:   PetscInt i;

449:   PetscFunctionBegin;
450:   PetscCheckSorted(n, X);
451:   *flg = PETSC_FALSE;
452:   for (i = 0; i < n - 1; i++) {
453:     if (X[i + 1] == X[i]) {
454:       *flg = PETSC_TRUE;
455:       break;
456:     }
457:   }
458:   PetscFunctionReturn(PETSC_SUCCESS);
459: }

461: /*@
462:   PetscSortedCheckDupsCount - Checks if a sorted `PetscCount` array has duplicates

464:   Not Collective

466:   Input Parameters:
467: + n - number of values
468: - X - sorted array of `PetscCount`

470:   Output Parameter:
471: . flg - True if the array has duplications, otherwise false

473:   Level: intermediate

475: .seealso: `PetscCount`, `PetscSortCount()`, `PetscSortedCheckDupsInt()`
476: @*/
477: PetscErrorCode PetscSortedCheckDupsCount(PetscCount n, const PetscCount X[], PetscBool *flg)
478: {
479:   PetscCount i;

481:   PetscFunctionBegin;
482:   PetscCheckSorted(n, X);
483:   *flg = PETSC_FALSE;
484:   for (i = 0; i < n - 1; i++) {
485:     if (X[i + 1] == X[i]) {
486:       *flg = PETSC_TRUE;
487:       break;
488:     }
489:   }
490:   PetscFunctionReturn(PETSC_SUCCESS);
491: }

493: /*@
494:   PetscSortRemoveDupsInt - Sorts an array of `PetscInt` in place in increasing order removes all duplicate entries

496:   Not Collective

498:   Input Parameters:
499: + n - number of values
500: - X - array of `PetscInt`

502:   Output Parameter:
503: . n - number of non-redundant values

505:   Level: intermediate

507: .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortedRemoveDupsInt()`
508: @*/
509: PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n, PetscInt X[])
510: {
511:   PetscFunctionBegin;
512:   PetscAssertPointer(n, 1);
513:   PetscCall(PetscSortInt(*n, X));
514:   PetscCall(PetscSortedRemoveDupsInt(n, X));
515:   PetscFunctionReturn(PETSC_SUCCESS);
516: }

518: /*@
519:   PetscFindInt - Finds the location of a `PetscInt` key in a sorted array of `PetscInt`

521:   Not Collective

523:   Input Parameters:
524: + key - the `PetscInt` key to locate
525: . n   - number of values in the array
526: - X   - array of `PetscInt`

528:   Output Parameter:
529: . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go

531:   Level: intermediate

533: .seealso: `PetscIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
534: @*/
535: PetscErrorCode PetscFindInt(PetscInt key, PetscCount n, const PetscInt X[], PetscInt *loc)
536: {
537:   PetscInt lo = 0, hi;

539:   PetscFunctionBegin;
540:   PetscAssertPointer(loc, 4);
541:   if (!n) {
542:     *loc = -1;
543:     PetscFunctionReturn(PETSC_SUCCESS);
544:   }
545:   PetscAssertPointer(X, 3);
546:   PetscCall(PetscIntCast(n, &hi));
547:   while (hi - lo > 1) {
548:     PetscInt mid = lo + (hi - lo) / 2;
549:     PetscAssert(X[lo] <= X[mid] && X[mid] <= X[hi - 1], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input array was not sorted: (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ")", X[lo], X[mid], X[hi - 1]);
550:     if (key < X[mid]) hi = mid;
551:     else lo = mid;
552:   }
553:   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
554:   PetscFunctionReturn(PETSC_SUCCESS);
555: }

557: /*@
558:   PetscFindCount - Finds the location of a `PetscCount` key in a sorted array of `PetscCount`

560:   Not Collective

562:   Input Parameters:
563: + key - the `PetscCount` key to locate
564: . n   - number of values in the array
565: - X   - array of `PetscCount`

567:   Output Parameter:
568: . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go

570:   Level: intermediate

572: .seealso: `PetscCount`, `PetscSortCount()`
573: @*/
574: PetscErrorCode PetscFindCount(PetscCount key, PetscCount n, const PetscCount X[], PetscCount *loc)
575: {
576:   PetscCount lo = 0, hi;

578:   PetscFunctionBegin;
579:   PetscAssertPointer(loc, 4);
580:   if (!n) {
581:     *loc = -1;
582:     PetscFunctionReturn(PETSC_SUCCESS);
583:   }
584:   PetscAssertPointer(X, 3);
585:   hi = n;
586:   while (hi - lo > 1) {
587:     PetscCount mid = lo + (hi - lo) / 2;
588:     PetscAssert(X[lo] <= X[mid] && X[mid] <= X[hi - 1], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input array was not sorted: (%" PetscCount_FMT ", %" PetscCount_FMT ", %" PetscCount_FMT ")", X[lo], X[mid], X[hi - 1]);
589:     if (key < X[mid]) hi = mid;
590:     else lo = mid;
591:   }
592:   *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
593:   PetscFunctionReturn(PETSC_SUCCESS);
594: }

596: /*@
597:   PetscCheckDupsInt - Checks if an `PetscInt` array has duplicates

599:   Not Collective

601:   Input Parameters:
602: + n - number of values in the array
603: - X - array of `PetscInt`

605:   Output Parameter:
606: . dups - True if the array has dups, otherwise false

608:   Level: intermediate

610: .seealso: `PetscSortRemoveDupsInt()`, `PetscSortedCheckDupsInt()`
611: @*/
612: PetscErrorCode PetscCheckDupsInt(PetscInt n, const PetscInt X[], PetscBool *dups)
613: {
614:   PetscHSetI ht;
615:   PetscBool  missing;

617:   PetscFunctionBegin;
618:   if (n) PetscAssertPointer(X, 2);
619:   PetscAssertPointer(dups, 3);
620:   *dups = PETSC_FALSE;
621:   if (n > 1) {
622:     PetscCall(PetscHSetICreate(&ht));
623:     PetscCall(PetscHSetIResize(ht, n));
624:     for (PetscInt i = 0; i < n; i++) {
625:       PetscCall(PetscHSetIQueryAdd(ht, X[i], &missing));
626:       if (!missing) {
627:         *dups = PETSC_TRUE;
628:         break;
629:       }
630:     }
631:     PetscCall(PetscHSetIDestroy(&ht));
632:   }
633:   PetscFunctionReturn(PETSC_SUCCESS);
634: }

636: /*@
637:   PetscFindMPIInt - Finds `PetscMPIInt` in a sorted array of `PetscMPIInt`

639:   Not Collective

641:   Input Parameters:
642: + key - the integer to locate
643: . n   - number of values in the array
644: - X   - array of `PetscMPIInt`

646:   Output Parameter:
647: . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go

649:   Level: intermediate

651: .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortInt()`, `PetscSortIntWithArray()`, `PetscSortRemoveDupsInt()`
652: @*/
653: PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscCount n, const PetscMPIInt X[], PetscInt *loc)
654: {
655:   PetscCount lo = 0, hi = n;

657:   PetscFunctionBegin;
658:   PetscAssertPointer(loc, 4);
659:   if (!n) {
660:     *loc = -1;
661:     PetscFunctionReturn(PETSC_SUCCESS);
662:   }
663:   PetscAssertPointer(X, 3);
664:   while (hi - lo > 1) {
665:     PetscCount mid = lo + (hi - lo) / 2;
666:     PetscAssert(X[lo] <= X[mid] && X[mid] <= X[hi - 1], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input array was not sorted: (%d, %d, %d)", X[lo], X[mid], X[hi - 1]);
667:     if (key < X[mid]) hi = mid;
668:     else lo = mid;
669:   }
670:   PetscCall(PetscIntCast(key == X[lo] ? lo : -(lo + (MPIU_Count)(key > X[lo]) + 1), loc));
671:   PetscFunctionReturn(PETSC_SUCCESS);
672: }

674: /*@
675:   PetscSortIntWithArray - Sorts an array of `PetscInt` in place in increasing order;
676:   changes a second array of `PetscInt` to match the sorted first array.

678:   Not Collective

680:   Input Parameters:
681: + n - number of values
682: . X - array of `PetscInt`
683: - Y - second array of `PetscInt`

685:   Level: intermediate

687: .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithCountArray()`
688: @*/
689: PetscErrorCode PetscSortIntWithArray(PetscCount n, PetscInt X[], PetscInt Y[])
690: {
691:   PetscInt pivot, t1, t2;

693:   PetscFunctionBegin;
694:   QuickSort2(PetscSortIntWithArray, X, Y, n, pivot, t1, t2);
695:   PetscFunctionReturn(PETSC_SUCCESS);
696: }

698: /*@
699:   PetscSortIntWithArrayPair - Sorts an array of `PetscInt` in place in increasing order;
700:   changes a pair of `PetscInt` arrays to match the sorted first array.

702:   Not Collective

704:   Input Parameters:
705: + n - number of values
706: . X - array of `PestcInt`
707: . Y - second array of `PestcInt` (first array of the pair)
708: - Z - third array of `PestcInt` (second array of the pair)

710:   Level: intermediate

712: .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithIntCountArrayPair()`
713: @*/
714: PetscErrorCode PetscSortIntWithArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscInt Z[])
715: {
716:   PetscInt pivot, t1, t2, t3;

718:   PetscFunctionBegin;
719:   QuickSort3(PetscSortIntWithArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
720:   PetscFunctionReturn(PETSC_SUCCESS);
721: }

723: /*@
724:   PetscSortIntWithMPIIntArray - Sorts an array of `PetscInt` in place in increasing order;
725:   changes a second array of `PetscMPI` to match the sorted first array.

727:   Not Collective

729:   Input Parameters:
730: + n - number of values
731: . X - array of `PetscInt`
732: - Y - second array of `PetscMPIInt`

734:   Level: intermediate

736: .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
737: @*/
738: PetscErrorCode PetscSortIntWithMPIIntArray(PetscCount n, PetscInt X[], PetscMPIInt Y[])
739: {
740:   PetscInt    pivot, t1;
741:   PetscMPIInt t2;

743:   PetscFunctionBegin;
744:   QuickSort2(PetscSortIntWithMPIIntArray, X, Y, n, pivot, t1, t2);
745:   PetscFunctionReturn(PETSC_SUCCESS);
746: }

748: /*@
749:   PetscSortIntWithCountArray - Sorts an array of `PetscInt` in place in increasing order;
750:   changes a second array of `PetscCount` to match the sorted first array.

752:   Not Collective

754:   Input Parameters:
755: + n - number of values
756: . X - array of `PetscInt`
757: - Y - second array of `PetscCount`

759:   Level: intermediate

761: .seealso: `PetscIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
762: @*/
763: PetscErrorCode PetscSortIntWithCountArray(PetscCount n, PetscInt X[], PetscCount Y[])
764: {
765:   PetscInt   pivot, t1;
766:   PetscCount t2;

768:   PetscFunctionBegin;
769:   QuickSort2(PetscSortIntWithCountArray, X, Y, n, pivot, t1, t2);
770:   PetscFunctionReturn(PETSC_SUCCESS);
771: }

773: /*@
774:   PetscSortIntWithIntCountArrayPair - Sorts an array of `PetscInt` in place in increasing order;
775:   changes a `PetscInt`  array and a `PetscCount` array to match the sorted first array.

777:   Not Collective

779:   Input Parameters:
780: + n - number of values
781: . X - array of `PetscInt`
782: . Y - second array of `PetscInt` (first array of the pair)
783: - Z - third array of `PetscCount` (second array of the pair)

785:   Level: intermediate

787:   Note:
788:   Usually X, Y are matrix row/column indices, and Z is a permutation array and therefore Z's type is PetscCount to allow 2B+ nonzeros even with 32-bit PetscInt.

790: .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortIntWithArray()`, `PetscIntSortSemiOrdered()`, `PetscSortIntWithArrayPair()`
791: @*/
792: PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount n, PetscInt X[], PetscInt Y[], PetscCount Z[])
793: {
794:   PetscInt   pivot, t1, t2; /* pivot is take from X[], so its type is still PetscInt */
795:   PetscCount t3;            /* temp for Z[] */

797:   PetscFunctionBegin;
798:   QuickSort3(PetscSortIntWithIntCountArrayPair, X, Y, Z, n, pivot, t1, t2, t3);
799:   PetscFunctionReturn(PETSC_SUCCESS);
800: }

802: /*@
803:   PetscSortedMPIInt - Determines whether the `PetscMPIInt` array is sorted.

805:   Not Collective

807:   Input Parameters:
808: + n - number of values
809: - X - array of `PetscMPIInt`

811:   Output Parameter:
812: . sorted - flag whether the array is sorted

814:   Level: intermediate

816: .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortMPIInt()`, `PetscSortedInt()`, `PetscSortedReal()`
817: @*/
818: PetscErrorCode PetscSortedMPIInt(PetscCount n, const PetscMPIInt X[], PetscBool *sorted)
819: {
820:   PetscFunctionBegin;
821:   PetscSorted(n, X, *sorted);
822:   PetscFunctionReturn(PETSC_SUCCESS);
823: }

825: /*@
826:   PetscSortMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order.

828:   Not Collective

830:   Input Parameters:
831: + n - number of values
832: - X - array of `PetscMPIInt`

834:   Level: intermediate

836:   Note:
837:   This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
838:   is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their
839:   code to see which routine is fastest.

841: .seealso: `PetscMPIIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`
842: @*/
843: PetscErrorCode PetscSortMPIInt(PetscCount n, PetscMPIInt X[])
844: {
845:   PetscMPIInt pivot, t1;

847:   PetscFunctionBegin;
848:   QuickSort1(PetscSortMPIInt, X, n, pivot, t1);
849:   PetscFunctionReturn(PETSC_SUCCESS);
850: }

852: /*@
853:   PetscSortRemoveDupsMPIInt - Sorts an array of `PetscMPIInt` in place in increasing order removes all duplicate entries

855:   Not Collective

857:   Input Parameters:
858: + n - number of values
859: - X - array of `PetscMPIInt`

861:   Output Parameter:
862: . n - number of non-redundant values

864:   Level: intermediate

866: .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
867: @*/
868: PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n, PetscMPIInt X[])
869: {
870:   PetscInt s = 0, N = *n, b = 0;

872:   PetscFunctionBegin;
873:   PetscCall(PetscSortMPIInt(N, X));
874:   for (PetscInt i = 0; i < N - 1; i++) {
875:     if (X[b + s + 1] != X[b]) {
876:       X[b + 1] = X[b + s + 1];
877:       b++;
878:     } else s++;
879:   }
880:   *n = N - s;
881:   PetscFunctionReturn(PETSC_SUCCESS);
882: }

884: /*@
885:   PetscSortMPIIntWithArray - Sorts an array of `PetscMPIInt` in place in increasing order;
886:   changes a second `PetscMPIInt` array to match the sorted first array.

888:   Not Collective

890:   Input Parameters:
891: + n - number of values
892: . X - array of `PetscMPIInt`
893: - Y - second array of `PetscMPIInt`

895:   Level: intermediate

897: .seealso: `PetscMPIIntSortSemiOrderedWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`
898: @*/
899: PetscErrorCode PetscSortMPIIntWithArray(PetscCount n, PetscMPIInt X[], PetscMPIInt Y[])
900: {
901:   PetscMPIInt pivot, t1, t2;

903:   PetscFunctionBegin;
904:   QuickSort2(PetscSortMPIIntWithArray, X, Y, n, pivot, t1, t2);
905:   PetscFunctionReturn(PETSC_SUCCESS);
906: }

908: /*@
909:   PetscSortMPIIntWithIntArray - Sorts an array of `PetscMPIInt` in place in increasing order;
910:   changes a second array of `PetscInt` to match the sorted first array.

912:   Not Collective

914:   Input Parameters:
915: + n - number of values
916: . X - array of `PetscMPIInt`
917: - Y - second array of `PetscInt`

919:   Level: intermediate

921:   Note:
922:   This routine is useful when one needs to sort MPI ranks with other integer arrays.

924: .seealso: `PetscSortMPIIntWithArray()`, `PetscIntSortSemiOrderedWithArray()`, `PetscTimSortWithArray()`
925: @*/
926: PetscErrorCode PetscSortMPIIntWithIntArray(PetscCount n, PetscMPIInt X[], PetscInt Y[])
927: {
928:   PetscMPIInt pivot, t1;
929:   PetscInt    t2;

931:   PetscFunctionBegin;
932:   QuickSort2(PetscSortMPIIntWithIntArray, X, Y, n, pivot, t1, t2);
933:   PetscFunctionReturn(PETSC_SUCCESS);
934: }

936: /*@
937:   PetscSortIntWithScalarArray - Sorts an array of `PetscInt` in place in increasing order;
938:   changes a second `PetscScalar` array to match the sorted first array.

940:   Not Collective

942:   Input Parameters:
943: + n - number of values
944: . X - array of `PetscInt`
945: - Y - second array of `PetscScalar`

947:   Level: intermediate

949: .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
950: @*/
951: PetscErrorCode PetscSortIntWithScalarArray(PetscCount n, PetscInt X[], PetscScalar Y[])
952: {
953:   PetscInt    pivot, t1;
954:   PetscScalar t2;

956:   PetscFunctionBegin;
957:   QuickSort2(PetscSortIntWithScalarArray, X, Y, n, pivot, t1, t2);
958:   PetscFunctionReturn(PETSC_SUCCESS);
959: }

961: /*@C
962:   PetscSortIntWithDataArray - Sorts an array of `PetscInt` in place in increasing order;
963:   changes a second array to match the sorted first INTEGER array.  Unlike other sort routines, the user must
964:   provide workspace (the size of an element in the data array) to use when sorting.

966:   Not Collective, No Fortran Support

968:   Input Parameters:
969: + n    - number of values
970: . X    - array of `PetscInt`
971: . Y    - second array of data
972: . size - sizeof elements in the data array in bytes
973: - t2   - workspace of "size" bytes used when sorting

975:   Level: intermediate

977: .seealso: `PetscTimSortWithArray()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
978: @*/
979: PetscErrorCode PetscSortIntWithDataArray(PetscCount n, PetscInt X[], void *Y, size_t size, void *t2)
980: {
981:   char      *YY = (char *)Y;
982:   PetscCount hi = n - 1;
983:   PetscInt   pivot, t1;

985:   PetscFunctionBegin;
986:   if (n < 8) {
987:     for (PetscCount i = 0; i < n; i++) {
988:       pivot = X[i];
989:       for (PetscCount j = i + 1; j < n; j++) {
990:         if (pivot > X[j]) {
991:           SWAP2Data(X[i], X[j], YY + size * i, YY + size * j, t1, t2, size);
992:           pivot = X[i];
993:         }
994:       }
995:     }
996:   } else {
997:     /* Two way partition */
998:     PetscCount l = 0, r = hi;

1000:     pivot = X[MEDIAN(X, hi)];
1001:     while (1) {
1002:       while (X[l] < pivot) l++;
1003:       while (X[r] > pivot) r--;
1004:       if (l >= r) {
1005:         r++;
1006:         break;
1007:       }
1008:       SWAP2Data(X[l], X[r], YY + size * l, YY + size * r, t1, t2, size);
1009:       l++;
1010:       r--;
1011:     }
1012:     PetscCall(PetscSortIntWithDataArray(l, X, Y, size, t2));
1013:     PetscCall(PetscSortIntWithDataArray(hi - r + 1, X + r, YY + size * r, size, t2));
1014:   }
1015:   PetscFunctionReturn(PETSC_SUCCESS);
1016: }

1018: /*@
1019:   PetscMergeIntArray -     Merges two SORTED `PetscInt` arrays, removes duplicate elements.

1021:   Not Collective

1023:   Input Parameters:
1024: + an - number of values in the first array
1025: . aI - first sorted array of `PetscInt`
1026: . bn - number of values in the second array
1027: - bI - second array of `PetscInt`

1029:   Output Parameters:
1030: + n - number of values in the merged array
1031: - L - merged sorted array, this is allocated if an array is not provided

1033:   Level: intermediate

1035: .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1036: @*/
1037: PetscErrorCode PetscMergeIntArray(PetscInt an, const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
1038: {
1039:   PetscInt *L_ = *L, ak, bk, k;

1041:   PetscFunctionBegin;
1042:   if (!L_) {
1043:     PetscCall(PetscMalloc1(an + bn, L));
1044:     L_ = *L;
1045:   }
1046:   k = ak = bk = 0;
1047:   while (ak < an && bk < bn) {
1048:     if (aI[ak] == bI[bk]) {
1049:       L_[k] = aI[ak];
1050:       ++ak;
1051:       ++bk;
1052:       ++k;
1053:     } else if (aI[ak] < bI[bk]) {
1054:       L_[k] = aI[ak];
1055:       ++ak;
1056:       ++k;
1057:     } else {
1058:       L_[k] = bI[bk];
1059:       ++bk;
1060:       ++k;
1061:     }
1062:   }
1063:   if (ak < an) {
1064:     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
1065:     k += (an - ak);
1066:   }
1067:   if (bk < bn) {
1068:     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
1069:     k += (bn - bk);
1070:   }
1071:   *n = k;
1072:   PetscFunctionReturn(PETSC_SUCCESS);
1073: }

1075: /*@
1076:   PetscMergeIntArrayPair -     Merges two SORTED `PetscInt` arrays that share NO common values along with an additional array of `PetscInt`.
1077:   The additional arrays are the same length as sorted arrays and are merged
1078:   in the order determined by the merging of the sorted pair.

1080:   Not Collective

1082:   Input Parameters:
1083: + an - number of values in the first array
1084: . aI - first sorted array of `PetscInt`
1085: . aJ - first additional array of `PetscInt`
1086: . bn - number of values in the second array
1087: . bI - second array of `PetscInt`
1088: - bJ - second additional of `PetscInt`

1090:   Output Parameters:
1091: + n - number of values in the merged array (== an + bn)
1092: . L - merged sorted array
1093: - J - merged additional array

1095:   Note:
1096:   if L or J point to non-null arrays then this routine will assume they are of the appropriate size and use them, otherwise this routine will allocate space for them

1098:   Level: intermediate

1100: .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1101: @*/
1102: PetscErrorCode PetscMergeIntArrayPair(PetscInt an, const PetscInt aI[], const PetscInt aJ[], PetscInt bn, const PetscInt bI[], const PetscInt bJ[], PetscInt *n, PetscInt *L[], PetscInt *J[])
1103: {
1104:   PetscInt n_, *L_, *J_, ak, bk, k;

1106:   PetscFunctionBegin;
1107:   PetscAssertPointer(L, 8);
1108:   PetscAssertPointer(J, 9);
1109:   n_ = an + bn;
1110:   *n = n_;
1111:   if (!*L) PetscCall(PetscMalloc1(n_, L));
1112:   L_ = *L;
1113:   if (!*J) PetscCall(PetscMalloc1(n_, J));
1114:   J_ = *J;
1115:   k = ak = bk = 0;
1116:   while (ak < an && bk < bn) {
1117:     if (aI[ak] <= bI[bk]) {
1118:       L_[k] = aI[ak];
1119:       J_[k] = aJ[ak];
1120:       ++ak;
1121:       ++k;
1122:     } else {
1123:       L_[k] = bI[bk];
1124:       J_[k] = bJ[bk];
1125:       ++bk;
1126:       ++k;
1127:     }
1128:   }
1129:   if (ak < an) {
1130:     PetscCall(PetscArraycpy(L_ + k, aI + ak, an - ak));
1131:     PetscCall(PetscArraycpy(J_ + k, aJ + ak, an - ak));
1132:     k += (an - ak);
1133:   }
1134:   if (bk < bn) {
1135:     PetscCall(PetscArraycpy(L_ + k, bI + bk, bn - bk));
1136:     PetscCall(PetscArraycpy(J_ + k, bJ + bk, bn - bk));
1137:   }
1138:   PetscFunctionReturn(PETSC_SUCCESS);
1139: }

1141: /*@
1142:   PetscMergeMPIIntArray -     Merges two SORTED `PetscMPIInt` arrays.

1144:   Not Collective

1146:   Input Parameters:
1147: + an - number of values in the first array
1148: . aI - first sorted array of `PetscMPIInt`
1149: . bn - number of values in the second array
1150: - bI - second array of `PetscMPIInt`

1152:   Output Parameters:
1153: + n - number of values in the merged array (<= an + bn)
1154: - L - merged sorted array, allocated if address of NULL pointer is passed

1156:   Level: intermediate

1158: .seealso: `PetscIntSortSemiOrdered()`, `PetscSortReal()`, `PetscSortIntWithPermutation()`, `PetscSortInt()`, `PetscSortIntWithArray()`
1159: @*/
1160: PetscErrorCode PetscMergeMPIIntArray(PetscInt an, const PetscMPIInt aI[], PetscInt bn, const PetscMPIInt bI[], PetscInt *n, PetscMPIInt **L)
1161: {
1162:   PetscInt ai, bi, k;

1164:   PetscFunctionBegin;
1165:   if (!*L) PetscCall(PetscMalloc1(an + bn, L));
1166:   for (ai = 0, bi = 0, k = 0; ai < an || bi < bn;) {
1167:     PetscInt t = -1;
1168:     for (; ai < an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
1169:     for (; bi < bn && bI[bi] == t; bi++);
1170:     for (; bi < bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
1171:     for (; ai < an && aI[ai] == t; ai++);
1172:   }
1173:   *n = k;
1174:   PetscFunctionReturn(PETSC_SUCCESS);
1175: }

1177: /*@C
1178:   PetscProcessTree - Prepares tree data to be displayed graphically

1180:   Not Collective, No Fortran Support

1182:   Input Parameters:
1183: + n        - number of values
1184: . mask     - indicates those entries in the tree, location 0 is always masked
1185: - parentid - indicates the parent of each entry

1187:   Output Parameters:
1188: + Nlevels   - the number of levels
1189: . Level     - for each node tells its level
1190: . Levelcnt  - the number of nodes on each level
1191: . Idbylevel - a list of ids on each of the levels, first level followed by second etc
1192: - Column    - for each id tells its column index

1194:   Level: developer

1196:   Note:
1197:   This code is not currently used

1199: .seealso: `PetscSortReal()`, `PetscSortIntWithPermutation()`
1200: @*/
1201: PetscErrorCode PetscProcessTree(PetscInt n, const PetscBool mask[], const PetscInt parentid[], PetscInt *Nlevels, PetscInt **Level, PetscInt **Levelcnt, PetscInt **Idbylevel, PetscInt **Column)
1202: {
1203:   PetscInt  i, j, cnt, nmask = 0, nlevels = 0, *level, *levelcnt, levelmax = 0, *workid, *workparentid, tcnt = 0, *idbylevel, *column;
1204:   PetscBool done = PETSC_FALSE;

1206:   PetscFunctionBegin;
1207:   PetscCheck(mask[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "Mask of 0th location must be set");
1208:   for (i = 0; i < n; i++) {
1209:     if (mask[i]) continue;
1210:     PetscCheck(parentid[i] != i, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Node labeled as own parent");
1211:     PetscCheck(!parentid[i] || !mask[parentid[i]], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Parent is masked");
1212:   }

1214:   for (i = 0; i < n; i++) {
1215:     if (!mask[i]) nmask++;
1216:   }

1218:   /* determine the level in the tree of each node */
1219:   PetscCall(PetscCalloc1(n, &level));

1221:   level[0] = 1;
1222:   while (!done) {
1223:     done = PETSC_TRUE;
1224:     for (i = 0; i < n; i++) {
1225:       if (mask[i]) continue;
1226:       if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
1227:       else if (!level[i]) done = PETSC_FALSE;
1228:     }
1229:   }
1230:   for (i = 0; i < n; i++) {
1231:     level[i]--;
1232:     nlevels = PetscMax(nlevels, level[i]);
1233:   }

1235:   /* count the number of nodes on each level and its max */
1236:   PetscCall(PetscCalloc1(nlevels, &levelcnt));
1237:   for (i = 0; i < n; i++) {
1238:     if (mask[i]) continue;
1239:     levelcnt[level[i] - 1]++;
1240:   }
1241:   for (i = 0; i < nlevels; i++) levelmax = PetscMax(levelmax, levelcnt[i]);

1243:   /* for each level sort the ids by the parent id */
1244:   PetscCall(PetscMalloc2(levelmax, &workid, levelmax, &workparentid));
1245:   PetscCall(PetscMalloc1(nmask, &idbylevel));
1246:   for (j = 1; j <= nlevels; j++) {
1247:     cnt = 0;
1248:     for (i = 0; i < n; i++) {
1249:       if (mask[i]) continue;
1250:       if (level[i] != j) continue;
1251:       workid[cnt]         = i;
1252:       workparentid[cnt++] = parentid[i];
1253:     }
1254:     /*  PetscIntView(cnt,workparentid,0);
1255:     PetscIntView(cnt,workid,0);
1256:     PetscCall(PetscSortIntWithArray(cnt,workparentid,workid));
1257:     PetscIntView(cnt,workparentid,0);
1258:     PetscIntView(cnt,workid,0);*/
1259:     PetscCall(PetscArraycpy(idbylevel + tcnt, workid, cnt));
1260:     tcnt += cnt;
1261:   }
1262:   PetscCheck(tcnt == nmask, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent count of unmasked nodes");
1263:   PetscCall(PetscFree2(workid, workparentid));

1265:   /* for each node list its column */
1266:   PetscCall(PetscMalloc1(n, &column));
1267:   cnt = 0;
1268:   for (j = 0; j < nlevels; j++) {
1269:     for (i = 0; i < levelcnt[j]; i++) column[idbylevel[cnt++]] = i;
1270:   }

1272:   *Nlevels   = nlevels;
1273:   *Level     = level;
1274:   *Levelcnt  = levelcnt;
1275:   *Idbylevel = idbylevel;
1276:   *Column    = column;
1277:   PetscFunctionReturn(PETSC_SUCCESS);
1278: }

1280: /*@
1281:   PetscParallelSortedInt - Check whether a `PetscInt` array, distributed over a communicator, is globally sorted.

1283:   Collective

1285:   Input Parameters:
1286: + comm - the MPI communicator
1287: . n    - the local number of `PetscInt`
1288: - keys - the local array of `PetscInt`

1290:   Output Parameters:
1291: . is_sorted - whether the array is globally sorted

1293:   Level: developer

1295: .seealso: `PetscParallelSortInt()`
1296: @*/
1297: PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
1298: {
1299:   PetscBool   sorted;
1300:   PetscInt    i, min, max, prevmax;
1301:   PetscMPIInt rank;

1303:   PetscFunctionBegin;
1304:   sorted = PETSC_TRUE;
1305:   min    = PETSC_INT_MAX;
1306:   max    = PETSC_INT_MIN;
1307:   if (n) {
1308:     min = keys[0];
1309:     max = keys[0];
1310:   }
1311:   for (i = 1; i < n; i++) {
1312:     if (keys[i] < keys[i - 1]) break;
1313:     min = PetscMin(min, keys[i]);
1314:     max = PetscMax(max, keys[i]);
1315:   }
1316:   if (i < n) sorted = PETSC_FALSE;
1317:   prevmax = PETSC_INT_MIN;
1318:   PetscCallMPI(MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm));
1319:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1320:   if (rank == 0) prevmax = PETSC_INT_MIN;
1321:   if (prevmax > min) sorted = PETSC_FALSE;
1322:   PetscCallMPI(MPIU_Allreduce(&sorted, is_sorted, 1, MPI_C_BOOL, MPI_LAND, comm));
1323:   PetscFunctionReturn(PETSC_SUCCESS);
1324: }