Actual source code: sorti.c
petsc-3.14.6 2021-03-30
2: /*
3: This file contains routines for sorting integers. Values are sorted in place.
4: One can use src/sys/tests/ex52.c for benchmarking.
5: */
6: #include <petsc/private/petscimpl.h>
7: #include <petsc/private/hashseti.h>
9: #define MEDIAN3(v,a,b,c) \
10: (v[a]<v[b] \
11: ? (v[b]<v[c] \
12: ? (b) \
13: : (v[a]<v[c] ? (c) : (a))) \
14: : (v[c]<v[b] \
15: ? (b) \
16: : (v[a]<v[c] ? (a) : (c))))
18: #define MEDIAN(v,right) MEDIAN3(v,right/4,right/2,right/4*3)
20: /* Swap one, two or three pairs. Each pair can have its own type */
21: #define SWAP1(a,b,t1) do {t1=a;a=b;b=t1;} while (0)
22: #define SWAP2(a,b,c,d,t1,t2) do {t1=a;a=b;b=t1; t2=c;c=d;d=t2;} while (0)
23: #define SWAP3(a,b,c,d,e,f,t1,t2,t3) do {t1=a;a=b;b=t1; t2=c;c=d;d=t2; t3=e;e=f;f=t3;} while (0)
25: /* Swap a & b, *c & *d. c, d, t2 are pointers to a type of size <siz> */
26: #define SWAP2Data(a,b,c,d,t1,t2,siz) \
27: do { \
29: t1=a; a=b; b=t1; \
30: PetscMemcpy(t2,c,siz); \
31: PetscMemcpy(c,d,siz); \
32: PetscMemcpy(d,t2,siz); \
33: } while (0)
35: /*
36: Partition X[lo,hi] into two parts: X[lo,l) <= pivot; X[r,hi] > pivot
38: Input Parameters:
39: + X - array to partition
40: . pivot - a pivot of X[]
41: . t1 - temp variable for X
42: - lo,hi - lower and upper bound of the array
44: Output Parameters:
45: + l,r - of type PetscInt
47: Notes:
48: The TwoWayPartition2/3 variants also partition other arrays along with X.
49: These arrays can have different types, so they provide their own temp t2,t3
50: */
51: #define TwoWayPartition1(X,pivot,t1,lo,hi,l,r) \
52: do { \
53: l = lo; \
54: r = hi; \
55: while (1) { \
56: while (X[l] < pivot) l++; \
57: while (X[r] > pivot) r--; \
58: if (l >= r) {r++; break;} \
59: SWAP1(X[l],X[r],t1); \
60: l++; \
61: r--; \
62: } \
63: } while (0)
65: /*
66: Partition X[lo,hi] into two parts: X[lo,l) >= pivot; X[r,hi] < pivot
68: Input Parameters:
69: + X - array to partition
70: . pivot - a pivot of X[]
71: . t1 - temp variable for X
72: - lo,hi - lower and upper bound of the array
74: Output Parameters:
75: + l,r - of type PetscInt
77: Notes:
78: The TwoWayPartition2/3 variants also partition other arrays along with X.
79: These arrays can have different types, so they provide their own temp t2,t3
80: */
81: #define TwoWayPartitionReverse1(X,pivot,t1,lo,hi,l,r) \
82: do { \
83: l = lo; \
84: r = hi; \
85: while (1) { \
86: while (X[l] > pivot) l++; \
87: while (X[r] < pivot) r--; \
88: if (l >= r) {r++; break;} \
89: SWAP1(X[l],X[r],t1); \
90: l++; \
91: r--; \
92: } \
93: } while (0)
95: #define TwoWayPartition2(X,Y,pivot,t1,t2,lo,hi,l,r) \
96: do { \
97: l = lo; \
98: r = hi; \
99: while (1) { \
100: while (X[l] < pivot) l++; \
101: while (X[r] > pivot) r--; \
102: if (l >= r) {r++; break;} \
103: SWAP2(X[l],X[r],Y[l],Y[r],t1,t2); \
104: l++; \
105: r--; \
106: } \
107: } while (0)
109: #define TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,lo,hi,l,r) \
110: do { \
111: l = lo; \
112: r = hi; \
113: while (1) { \
114: while (X[l] < pivot) l++; \
115: while (X[r] > pivot) r--; \
116: if (l >= r) {r++; break;} \
117: SWAP3(X[l],X[r],Y[l],Y[r],Z[l],Z[r],t1,t2,t3); \
118: l++; \
119: r--; \
120: } \
121: } while (0)
123: /* Templates for similar functions used below */
124: #define QuickSort1(FuncName,X,n,pivot,t1,ierr) \
125: do { \
126: PetscInt i,j,p,l,r,hi=n-1; \
127: if (n < 8) { \
128: for (i=0; i<n; i++) { \
129: pivot = X[i]; \
130: for (j=i+1; j<n; j++) { \
131: if (pivot > X[j]) { \
132: SWAP1(X[i],X[j],t1); \
133: pivot = X[i]; \
134: } \
135: } \
136: } \
137: } else { \
138: p = MEDIAN(X,hi); \
139: pivot = X[p]; \
140: TwoWayPartition1(X,pivot,t1,0,hi,l,r); \
141: FuncName(l,X); \
142: FuncName(hi-r+1,X+r); \
143: } \
144: } while (0)
146: /* Templates for similar functions used below */
147: #define QuickSortReverse1(FuncName,X,n,pivot,t1,ierr) \
148: do { \
149: PetscInt i,j,p,l,r,hi=n-1; \
150: if (n < 8) { \
151: for (i=0; i<n; i++) { \
152: pivot = X[i]; \
153: for (j=i+1; j<n; j++) { \
154: if (pivot < X[j]) { \
155: SWAP1(X[i],X[j],t1); \
156: pivot = X[i]; \
157: } \
158: } \
159: } \
160: } else { \
161: p = MEDIAN(X,hi); \
162: pivot = X[p]; \
163: TwoWayPartitionReverse1(X,pivot,t1,0,hi,l,r); \
164: FuncName(l,X); \
165: FuncName(hi-r+1,X+r); \
166: } \
167: } while (0)
169: #define QuickSort2(FuncName,X,Y,n,pivot,t1,t2,ierr) \
170: do { \
171: PetscInt i,j,p,l,r,hi=n-1; \
172: if (n < 8) { \
173: for (i=0; i<n; i++) { \
174: pivot = X[i]; \
175: for (j=i+1; j<n; j++) { \
176: if (pivot > X[j]) { \
177: SWAP2(X[i],X[j],Y[i],Y[j],t1,t2); \
178: pivot = X[i]; \
179: } \
180: } \
181: } \
182: } else { \
183: p = MEDIAN(X,hi); \
184: pivot = X[p]; \
185: TwoWayPartition2(X,Y,pivot,t1,t2,0,hi,l,r); \
186: FuncName(l,X,Y); \
187: FuncName(hi-r+1,X+r,Y+r); \
188: } \
189: } while (0)
191: #define QuickSort3(FuncName,X,Y,Z,n,pivot,t1,t2,t3,ierr) \
192: do { \
193: PetscInt i,j,p,l,r,hi=n-1; \
194: if (n < 8) { \
195: for (i=0; i<n; i++) { \
196: pivot = X[i]; \
197: for (j=i+1; j<n; j++) { \
198: if (pivot > X[j]) { \
199: SWAP3(X[i],X[j],Y[i],Y[j],Z[i],Z[j],t1,t2,t3); \
200: pivot = X[i]; \
201: } \
202: } \
203: } \
204: } else { \
205: p = MEDIAN(X,hi); \
206: pivot = X[p]; \
207: TwoWayPartition3(X,Y,Z,pivot,t1,t2,t3,0,hi,l,r); \
208: FuncName(l,X,Y,Z); \
209: FuncName(hi-r+1,X+r,Y+r,Z+r); \
210: } \
211: } while (0)
213: /*@
214: PetscSortedInt - Determines whether the array is sorted.
216: Not Collective
218: Input Parameters:
219: + n - number of values
220: - X - array of integers
222: Output Parameters:
223: . sorted - flag whether the array is sorted
225: Level: intermediate
227: .seealso: PetscSortInt(), PetscSortedMPIInt(), PetscSortedReal()
228: @*/
229: PetscErrorCode PetscSortedInt(PetscInt n,const PetscInt X[],PetscBool *sorted)
230: {
232: PetscSorted(n,X,*sorted);
233: return(0);
234: }
236: /*@
237: PetscSortInt - Sorts an array of integers in place in increasing order.
239: Not Collective
241: Input Parameters:
242: + n - number of values
243: - X - array of integers
245: Notes:
246: This function serves as an alternative to PetscIntSortSemiOrdered(), and may perform faster especially if the array
247: is completely random. There are exceptions to this and so it is __highly__ recomended that the user benchmark their
248: code to see which routine is fastest.
250: Level: intermediate
252: .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation()
253: @*/
254: PetscErrorCode PetscSortInt(PetscInt n,PetscInt X[])
255: {
257: PetscInt pivot,t1;
260: QuickSort1(PetscSortInt,X,n,pivot,t1,ierr);
261: return(0);
262: }
264: /*@
265: PetscSortReverseInt - Sorts an array of integers in place in decreasing order.
267: Not Collective
269: Input Parameters:
270: + n - number of values
271: - X - array of integers
273: Level: intermediate
275: .seealso: PetscIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithPermutation()
276: @*/
277: PetscErrorCode PetscSortReverseInt(PetscInt n,PetscInt X[])
278: {
280: PetscInt pivot,t1;
283: QuickSortReverse1(PetscSortReverseInt,X,n,pivot,t1,ierr);
284: return(0);
285: }
287: /*@
288: PetscSortedRemoveDupsInt - Removes all duplicate entries of a sorted input array
290: Not Collective
292: Input Parameters:
293: + n - number of values
294: - X - sorted array of integers
296: Output Parameter:
297: . n - number of non-redundant values
299: Level: intermediate
301: .seealso: PetscSortInt()
302: @*/
303: PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *n,PetscInt X[])
304: {
305: PetscInt i,s = 0,N = *n, b = 0;
309: for (i=0; i<N-1; i++) {
310: if (X[b+s+1] != X[b]) {
311: X[b+1] = X[b+s+1]; b++;
312: } else s++;
313: }
314: *n = N - s;
315: return(0);
316: }
318: /*@
319: PetscSortRemoveDupsInt - Sorts an array of integers in place in increasing order removes all duplicate entries
321: Not Collective
323: Input Parameters:
324: + n - number of values
325: - X - array of integers
327: Output Parameter:
328: . n - number of non-redundant values
330: Level: intermediate
332: .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt(), PetscSortedRemoveDupsInt()
333: @*/
334: PetscErrorCode PetscSortRemoveDupsInt(PetscInt *n,PetscInt X[])
335: {
339: PetscSortInt(*n,X);
340: PetscSortedRemoveDupsInt(n,X);
341: return(0);
342: }
344: /*@
345: PetscFindInt - Finds integer in a sorted array of integers
347: Not Collective
349: Input Parameters:
350: + key - the integer to locate
351: . n - number of values in the array
352: - X - array of integers
354: Output Parameter:
355: . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
357: Level: intermediate
359: .seealso: PetscIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
360: @*/
361: PetscErrorCode PetscFindInt(PetscInt key, PetscInt n, const PetscInt X[], PetscInt *loc)
362: {
363: PetscInt lo = 0,hi = n;
367: if (!n) {*loc = -1; return(0);}
370: while (hi - lo > 1) {
371: PetscInt mid = lo + (hi - lo)/2;
372: if (key < X[mid]) hi = mid;
373: else lo = mid;
374: }
375: *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
376: return(0);
377: }
379: /*@
382: Not Collective
384: Input Parameters:
385: + n - number of values in the array
386: - X - array of integers
389: Output Parameter:
390: . dups - True if the array has dups, otherwise false
392: Level: intermediate
394: .seealso: PetscSortRemoveDupsInt()
395: @*/
397: {
399: PetscInt i;
400: PetscHSetI ht;
401: PetscBool missing;
405: *dups = PETSC_FALSE;
406: if (n > 1) {
407: PetscHSetICreate(&ht);
408: PetscHSetIResize(ht,n);
409: for (i=0; i<n; i++) {
410: PetscHSetIQueryAdd(ht,X[i],&missing);
411: if (!missing) {*dups = PETSC_TRUE; break;}
412: }
413: PetscHSetIDestroy(&ht);
414: }
415: return(0);
416: }
418: /*@
419: PetscFindMPIInt - Finds MPI integer in a sorted array of integers
421: Not Collective
423: Input Parameters:
424: + key - the integer to locate
425: . n - number of values in the array
426: - X - array of integers
428: Output Parameter:
429: . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go
431: Level: intermediate
433: .seealso: PetscMPIIntSortSemiOrdered(), PetscSortInt(), PetscSortIntWithArray(), PetscSortRemoveDupsInt()
434: @*/
435: PetscErrorCode PetscFindMPIInt(PetscMPIInt key, PetscInt n, const PetscMPIInt X[], PetscInt *loc)
436: {
437: PetscInt lo = 0,hi = n;
441: if (!n) {*loc = -1; return(0);}
444: while (hi - lo > 1) {
445: PetscInt mid = lo + (hi - lo)/2;
446: if (key < X[mid]) hi = mid;
447: else lo = mid;
448: }
449: *loc = key == X[lo] ? lo : -(lo + (key > X[lo]) + 1);
450: return(0);
451: }
453: /*@
454: PetscSortIntWithArray - Sorts an array of integers in place in increasing order;
455: changes a second array to match the sorted first array.
457: Not Collective
459: Input Parameters:
460: + n - number of values
461: . X - array of integers
462: - Y - second array of integers
464: Level: intermediate
466: .seealso: PetscIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
467: @*/
468: PetscErrorCode PetscSortIntWithArray(PetscInt n,PetscInt X[],PetscInt Y[])
469: {
471: PetscInt pivot,t1,t2;
474: QuickSort2(PetscSortIntWithArray,X,Y,n,pivot,t1,t2,ierr);
475: return(0);
476: }
478: /*@
479: PetscSortIntWithArrayPair - Sorts an array of integers in place in increasing order;
480: changes a pair of integer arrays to match the sorted first array.
482: Not Collective
484: Input Parameters:
485: + n - number of values
486: . X - array of integers
487: . Y - second array of integers (first array of the pair)
488: - Z - third array of integers (second array of the pair)
490: Level: intermediate
492: .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortIntWithArray(), PetscIntSortSemiOrdered()
493: @*/
494: PetscErrorCode PetscSortIntWithArrayPair(PetscInt n,PetscInt X[],PetscInt Y[],PetscInt Z[])
495: {
497: PetscInt pivot,t1,t2,t3;
500: QuickSort3(PetscSortIntWithArrayPair,X,Y,Z,n,pivot,t1,t2,t3,ierr);
501: return(0);
502: }
504: /*@
505: PetscSortedMPIInt - Determines whether the array is sorted.
507: Not Collective
509: Input Parameters:
510: + n - number of values
511: - X - array of integers
513: Output Parameters:
514: . sorted - flag whether the array is sorted
516: Level: intermediate
518: .seealso: PetscMPIIntSortSemiOrdered(), PetscSortMPIInt(), PetscSortedInt(), PetscSortedReal()
519: @*/
520: PetscErrorCode PetscSortedMPIInt(PetscInt n,const PetscMPIInt X[],PetscBool *sorted)
521: {
523: PetscSorted(n,X,*sorted);
524: return(0);
525: }
527: /*@
528: PetscSortMPIInt - Sorts an array of MPI integers in place in increasing order.
530: Not Collective
532: Input Parameters:
533: + n - number of values
534: - X - array of integers
536: Level: intermediate
538: Notes:
539: This function serves as an alternative to PetscMPIIntSortSemiOrdered(), and may perform faster especially if the array
540: is completely random. There are exceptions to this and so it is __highly__ recomended that the user benchmark their
541: code to see which routine is fastest.
543: .seealso: PetscMPIIntSortSemiOrdered(), PetscSortReal(), PetscSortIntWithPermutation()
544: @*/
545: PetscErrorCode PetscSortMPIInt(PetscInt n,PetscMPIInt X[])
546: {
548: PetscMPIInt pivot,t1;
551: QuickSort1(PetscSortMPIInt,X,n,pivot,t1,ierr);
552: return(0);
553: }
555: /*@
556: PetscSortRemoveDupsMPIInt - Sorts an array of MPI integers in place in increasing order removes all duplicate entries
558: Not Collective
560: Input Parameters:
561: + n - number of values
562: - X - array of integers
564: Output Parameter:
565: . n - number of non-redundant values
567: Level: intermediate
569: .seealso: PetscSortReal(), PetscSortIntWithPermutation(), PetscSortInt()
570: @*/
571: PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *n,PetscMPIInt X[])
572: {
574: PetscInt i,s = 0,N = *n, b = 0;
577: PetscSortMPIInt(N,X);
578: for (i=0; i<N-1; i++) {
579: if (X[b+s+1] != X[b]) {
580: X[b+1] = X[b+s+1]; b++;
581: } else s++;
582: }
583: *n = N - s;
584: return(0);
585: }
587: /*@
588: PetscSortMPIIntWithArray - Sorts an array of integers in place in increasing order;
589: changes a second array to match the sorted first array.
591: Not Collective
593: Input Parameters:
594: + n - number of values
595: . X - array of integers
596: - Y - second array of integers
598: Level: intermediate
600: .seealso: PetscMPIIntSortSemiOrderedWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt()
601: @*/
602: PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt n,PetscMPIInt X[],PetscMPIInt Y[])
603: {
605: PetscMPIInt pivot,t1,t2;
608: QuickSort2(PetscSortMPIIntWithArray,X,Y,n,pivot,t1,t2,ierr);
609: return(0);
610: }
612: /*@
613: PetscSortMPIIntWithIntArray - Sorts an array of MPI integers in place in increasing order;
614: changes a second array of Petsc intergers to match the sorted first array.
616: Not Collective
618: Input Parameters:
619: + n - number of values
620: . X - array of MPI integers
621: - Y - second array of Petsc integers
623: Level: intermediate
625: Notes: this routine is useful when one needs to sort MPI ranks with other integer arrays.
627: .seealso: PetscSortMPIIntWithArray(), PetscIntSortSemiOrderedWithArray(), PetscTimSortWithArray()
628: @*/
629: PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt n,PetscMPIInt X[],PetscInt Y[])
630: {
632: PetscMPIInt pivot,t1;
633: PetscInt t2;
636: QuickSort2(PetscSortMPIIntWithIntArray,X,Y,n,pivot,t1,t2,ierr);
637: return(0);
638: }
640: /*@
641: PetscSortIntWithScalarArray - Sorts an array of integers in place in increasing order;
642: changes a second SCALAR array to match the sorted first INTEGER array.
644: Not Collective
646: Input Parameters:
647: + n - number of values
648: . X - array of integers
649: - Y - second array of scalars
651: Level: intermediate
653: .seealso: PetscTimSortWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
654: @*/
655: PetscErrorCode PetscSortIntWithScalarArray(PetscInt n,PetscInt X[],PetscScalar Y[])
656: {
658: PetscInt pivot,t1;
659: PetscScalar t2;
662: QuickSort2(PetscSortIntWithScalarArray,X,Y,n,pivot,t1,t2,ierr);
663: return(0);
664: }
666: /*@C
667: PetscSortIntWithDataArray - Sorts an array of integers in place in increasing order;
668: changes a second array to match the sorted first INTEGER array. Unlike other sort routines, the user must
669: provide workspace (the size of an element in the data array) to use when sorting.
671: Not Collective
673: Input Parameters:
674: + n - number of values
675: . X - array of integers
676: . Y - second array of data
677: . size - sizeof elements in the data array in bytes
678: - t2 - workspace of "size" bytes used when sorting
680: Level: intermediate
682: .seealso: PetscTimSortWithArray(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
683: @*/
684: PetscErrorCode PetscSortIntWithDataArray(PetscInt n,PetscInt X[],void *Y,size_t size,void *t2)
685: {
687: char *YY = (char*)Y;
688: PetscInt i,j,p,t1,pivot,hi=n-1,l,r;
691: if (n<8) {
692: for (i=0; i<n; i++) {
693: pivot = X[i];
694: for (j=i+1; j<n; j++) {
695: if (pivot > X[j]) {
696: SWAP2Data(X[i],X[j],YY+size*i,YY+size*j,t1,t2,size);
697: pivot = X[i];
698: }
699: }
700: }
701: } else {
702: /* Two way partition */
703: p = MEDIAN(X,hi);
704: pivot = X[p];
705: l = 0;
706: r = hi;
707: while (1) {
708: while (X[l] < pivot) l++;
709: while (X[r] > pivot) r--;
710: if (l >= r) {r++; break;}
711: SWAP2Data(X[l],X[r],YY+size*l,YY+size*r,t1,t2,size);
712: l++;
713: r--;
714: }
715: PetscSortIntWithDataArray(l,X,Y,size,t2);
716: PetscSortIntWithDataArray(hi-r+1,X+r,YY+size*r,size,t2);
717: }
718: return(0);
719: }
721: /*@
722: PetscMergeIntArray - Merges two SORTED integer arrays, removes duplicate elements.
724: Not Collective
726: Input Parameters:
727: + an - number of values in the first array
728: . aI - first sorted array of integers
729: . bn - number of values in the second array
730: - bI - second array of integers
732: Output Parameters:
733: + n - number of values in the merged array
734: - L - merged sorted array, this is allocated if an array is not provided
736: Level: intermediate
738: .seealso: PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
739: @*/
740: PetscErrorCode PetscMergeIntArray(PetscInt an,const PetscInt aI[], PetscInt bn, const PetscInt bI[], PetscInt *n, PetscInt **L)
741: {
743: PetscInt *L_ = *L, ak, bk, k;
745: if (!L_) {
746: PetscMalloc1(an+bn, L);
747: L_ = *L;
748: }
749: k = ak = bk = 0;
750: while (ak < an && bk < bn) {
751: if (aI[ak] == bI[bk]) {
752: L_[k] = aI[ak];
753: ++ak;
754: ++bk;
755: ++k;
756: } else if (aI[ak] < bI[bk]) {
757: L_[k] = aI[ak];
758: ++ak;
759: ++k;
760: } else {
761: L_[k] = bI[bk];
762: ++bk;
763: ++k;
764: }
765: }
766: if (ak < an) {
767: PetscArraycpy(L_+k,aI+ak,an-ak);
768: k += (an-ak);
769: }
770: if (bk < bn) {
771: PetscArraycpy(L_+k,bI+bk,bn-bk);
772: k += (bn-bk);
773: }
774: *n = k;
775: return(0);
776: }
778: /*@
779: PetscMergeIntArrayPair - Merges two SORTED integer arrays that share NO common values along with an additional array of integers.
780: The additional arrays are the same length as sorted arrays and are merged
781: in the order determined by the merging of the sorted pair.
783: Not Collective
785: Input Parameters:
786: + an - number of values in the first array
787: . aI - first sorted array of integers
788: . aJ - first additional array of integers
789: . bn - number of values in the second array
790: . bI - second array of integers
791: - bJ - second additional of integers
793: Output Parameters:
794: + n - number of values in the merged array (== an + bn)
795: . L - merged sorted array
796: - J - merged additional array
798: Notes:
799: if L or J point to non-null arrays then this routine will assume they are of the approproate size and use them, otherwise this routine will allocate space for them
800: Level: intermediate
802: .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
803: @*/
804: PetscErrorCode PetscMergeIntArrayPair(PetscInt an,const PetscInt aI[], const PetscInt aJ[], PetscInt bn, const PetscInt bI[], const PetscInt bJ[], PetscInt *n, PetscInt **L, PetscInt **J)
805: {
807: PetscInt n_, *L_, *J_, ak, bk, k;
812: n_ = an + bn;
813: *n = n_;
814: if (!*L) {
815: PetscMalloc1(n_, L);
816: }
817: L_ = *L;
818: if (!*J) {
819: PetscMalloc1(n_, J);
820: }
821: J_ = *J;
822: k = ak = bk = 0;
823: while (ak < an && bk < bn) {
824: if (aI[ak] <= bI[bk]) {
825: L_[k] = aI[ak];
826: J_[k] = aJ[ak];
827: ++ak;
828: ++k;
829: } else {
830: L_[k] = bI[bk];
831: J_[k] = bJ[bk];
832: ++bk;
833: ++k;
834: }
835: }
836: if (ak < an) {
837: PetscArraycpy(L_+k,aI+ak,an-ak);
838: PetscArraycpy(J_+k,aJ+ak,an-ak);
839: k += (an-ak);
840: }
841: if (bk < bn) {
842: PetscArraycpy(L_+k,bI+bk,bn-bk);
843: PetscArraycpy(J_+k,bJ+bk,bn-bk);
844: }
845: return(0);
846: }
848: /*@
849: PetscMergeMPIIntArray - Merges two SORTED integer arrays.
851: Not Collective
853: Input Parameters:
854: + an - number of values in the first array
855: . aI - first sorted array of integers
856: . bn - number of values in the second array
857: - bI - second array of integers
859: Output Parameters:
860: + n - number of values in the merged array (<= an + bn)
861: - L - merged sorted array, allocated if address of NULL pointer is passed
863: Level: intermediate
865: .seealso: PetscIntSortSemiOrdered(), PetscSortReal(), PetscSortIntPermutation(), PetscSortInt(), PetscSortIntWithArray()
866: @*/
867: PetscErrorCode PetscMergeMPIIntArray(PetscInt an,const PetscMPIInt aI[],PetscInt bn,const PetscMPIInt bI[],PetscInt *n,PetscMPIInt **L)
868: {
870: PetscInt ai,bi,k;
873: if (!*L) {PetscMalloc1((an+bn),L);}
874: for (ai=0,bi=0,k=0; ai<an || bi<bn;) {
875: PetscInt t = -1;
876: for (; ai<an && (!bn || aI[ai] <= bI[bi]); ai++) (*L)[k++] = t = aI[ai];
877: for (; bi<bn && bI[bi] == t; bi++);
878: for (; bi<bn && (!an || bI[bi] <= aI[ai]); bi++) (*L)[k++] = t = bI[bi];
879: for (; ai<an && aI[ai] == t; ai++);
880: }
881: *n = k;
882: return(0);
883: }
885: /*@C
886: PetscProcessTree - Prepares tree data to be displayed graphically
888: Not Collective
890: Input Parameters:
891: + n - number of values
892: . mask - indicates those entries in the tree, location 0 is always masked
893: - parentid - indicates the parent of each entry
895: Output Parameters:
896: + Nlevels - the number of levels
897: . Level - for each node tells its level
898: . Levelcnts - the number of nodes on each level
899: . Idbylevel - a list of ids on each of the levels, first level followed by second etc
900: - Column - for each id tells its column index
902: Level: developer
904: Notes:
905: This code is not currently used
907: .seealso: PetscSortReal(), PetscSortIntWithPermutation()
908: @*/
909: PetscErrorCode PetscProcessTree(PetscInt n,const PetscBool mask[],const PetscInt parentid[],PetscInt *Nlevels,PetscInt **Level,PetscInt **Levelcnt,PetscInt **Idbylevel,PetscInt **Column)
910: {
911: PetscInt i,j,cnt,nmask = 0,nlevels = 0,*level,*levelcnt,levelmax = 0,*workid,*workparentid,tcnt = 0,*idbylevel,*column;
913: PetscBool done = PETSC_FALSE;
916: if (!mask[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mask of 0th location must be set");
917: for (i=0; i<n; i++) {
918: if (mask[i]) continue;
919: if (parentid[i] == i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Node labeled as own parent");
920: if (parentid[i] && mask[parentid[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parent is masked");
921: }
923: for (i=0; i<n; i++) {
924: if (!mask[i]) nmask++;
925: }
927: /* determine the level in the tree of each node */
928: PetscCalloc1(n,&level);
930: level[0] = 1;
931: while (!done) {
932: done = PETSC_TRUE;
933: for (i=0; i<n; i++) {
934: if (mask[i]) continue;
935: if (!level[i] && level[parentid[i]]) level[i] = level[parentid[i]] + 1;
936: else if (!level[i]) done = PETSC_FALSE;
937: }
938: }
939: for (i=0; i<n; i++) {
940: level[i]--;
941: nlevels = PetscMax(nlevels,level[i]);
942: }
944: /* count the number of nodes on each level and its max */
945: PetscCalloc1(nlevels,&levelcnt);
946: for (i=0; i<n; i++) {
947: if (mask[i]) continue;
948: levelcnt[level[i]-1]++;
949: }
950: for (i=0; i<nlevels;i++) levelmax = PetscMax(levelmax,levelcnt[i]);
952: /* for each level sort the ids by the parent id */
953: PetscMalloc2(levelmax,&workid,levelmax,&workparentid);
954: PetscMalloc1(nmask,&idbylevel);
955: for (j=1; j<=nlevels;j++) {
956: cnt = 0;
957: for (i=0; i<n; i++) {
958: if (mask[i]) continue;
959: if (level[i] != j) continue;
960: workid[cnt] = i;
961: workparentid[cnt++] = parentid[i];
962: }
963: /* PetscIntView(cnt,workparentid,0);
964: PetscIntView(cnt,workid,0);
965: PetscSortIntWithArray(cnt,workparentid,workid);
966: PetscIntView(cnt,workparentid,0);
967: PetscIntView(cnt,workid,0);*/
968: PetscArraycpy(idbylevel+tcnt,workid,cnt);
969: tcnt += cnt;
970: }
971: if (tcnt != nmask) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent count of unmasked nodes");
972: PetscFree2(workid,workparentid);
974: /* for each node list its column */
975: PetscMalloc1(n,&column);
976: cnt = 0;
977: for (j=0; j<nlevels; j++) {
978: for (i=0; i<levelcnt[j]; i++) {
979: column[idbylevel[cnt++]] = i;
980: }
981: }
983: *Nlevels = nlevels;
984: *Level = level;
985: *Levelcnt = levelcnt;
986: *Idbylevel = idbylevel;
987: *Column = column;
988: return(0);
989: }
991: /*@
992: PetscParallelSortedInt - Check whether an integer array, distributed over a communicator, is globally sorted.
994: Collective
996: Input Parameters:
997: + comm - the MPI communicator
998: . n - the local number of integers
999: - keys - the local array of integers
1001: Output Parameters:
1002: . is_sorted - whether the array is globally sorted
1004: Level: developer
1006: .seealso: PetscParallelSortInt()
1007: @*/
1008: PetscErrorCode PetscParallelSortedInt(MPI_Comm comm, PetscInt n, const PetscInt keys[], PetscBool *is_sorted)
1009: {
1010: PetscBool sorted;
1011: PetscInt i, min, max, prevmax;
1012: PetscMPIInt rank;
1016: sorted = PETSC_TRUE;
1017: min = PETSC_MAX_INT;
1018: max = PETSC_MIN_INT;
1019: if (n) {
1020: min = keys[0];
1021: max = keys[0];
1022: }
1023: for (i = 1; i < n; i++) {
1024: if (keys[i] < keys[i - 1]) break;
1025: min = PetscMin(min,keys[i]);
1026: max = PetscMax(max,keys[i]);
1027: }
1028: if (i < n) sorted = PETSC_FALSE;
1029: prevmax = PETSC_MIN_INT;
1030: MPI_Exscan(&max, &prevmax, 1, MPIU_INT, MPI_MAX, comm);
1031: MPI_Comm_rank(comm, &rank);
1032: if (!rank) prevmax = PETSC_MIN_INT;
1033: if (prevmax > min) sorted = PETSC_FALSE;
1034: MPI_Allreduce(&sorted, is_sorted, 1, MPIU_BOOL, MPI_LAND, comm);
1035: return(0);
1036: }