xref: /petsc/src/vec/is/utils/psort.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1ce605777SToby Isaac 
2ce605777SToby Isaac #include <petsc/private/petscimpl.h>
3579e05a6SToby Isaac #include <petscis.h> /*I "petscis.h" I*/
4ce605777SToby Isaac 
5ce605777SToby Isaac /* This is the bitonic merge that works on non-power-of-2 sizes found at http://www.iti.fh-flensburg.de/lang/algorithmen/sortieren/bitonic/oddn.htm */
69371c9d4SSatish Balay static PetscErrorCode PetscParallelSortInt_Bitonic_Merge(MPI_Comm comm, PetscMPIInt tag, PetscMPIInt rankStart, PetscMPIInt rankEnd, PetscMPIInt rank, PetscMPIInt n, PetscInt keys[], PetscInt buffer[], PetscBool forward) {
7ce605777SToby Isaac   PetscInt diff;
8ce605777SToby Isaac   PetscInt split, mid, partner;
9ce605777SToby Isaac 
10ce605777SToby Isaac   PetscFunctionBegin;
11ce605777SToby Isaac   diff = rankEnd - rankStart;
12ce605777SToby Isaac   if (diff <= 0) PetscFunctionReturn(0);
13ce605777SToby Isaac   if (diff == 1) {
14ce605777SToby Isaac     if (forward) {
159566063dSJacob Faibussowitsch       PetscCall(PetscSortInt((PetscInt)n, keys));
16ce605777SToby Isaac     } else {
179566063dSJacob Faibussowitsch       PetscCall(PetscSortReverseInt((PetscInt)n, keys));
18ce605777SToby Isaac     }
19ce605777SToby Isaac     PetscFunctionReturn(0);
20ce605777SToby Isaac   }
21ce605777SToby Isaac   split = 1;
22ce605777SToby Isaac   while (2 * split < diff) split *= 2;
23ce605777SToby Isaac   mid = rankStart + split;
24ce605777SToby Isaac   if (rank < mid) {
25ce605777SToby Isaac     partner = rank + split;
26ce605777SToby Isaac   } else {
27ce605777SToby Isaac     partner = rank - split;
28ce605777SToby Isaac   }
29ce605777SToby Isaac   if (partner < rankEnd) {
30ce605777SToby Isaac     PetscMPIInt i;
31ce605777SToby Isaac 
329566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Sendrecv(keys, n, MPIU_INT, partner, tag, buffer, n, MPIU_INT, partner, tag, comm, MPI_STATUS_IGNORE));
33ce605777SToby Isaac     if ((rank < partner) == (forward == PETSC_TRUE)) {
349371c9d4SSatish Balay       for (i = 0; i < n; i++) { keys[i] = (keys[i] <= buffer[i]) ? keys[i] : buffer[i]; }
35ce605777SToby Isaac     } else {
369371c9d4SSatish Balay       for (i = 0; i < n; i++) { keys[i] = (keys[i] > buffer[i]) ? keys[i] : buffer[i]; }
37ce605777SToby Isaac     }
38ce605777SToby Isaac   }
39ce605777SToby Isaac   /* divide and conquer */
40ce605777SToby Isaac   if (rank < mid) {
419566063dSJacob Faibussowitsch     PetscCall(PetscParallelSortInt_Bitonic_Merge(comm, tag, rankStart, mid, rank, n, keys, buffer, forward));
42ce605777SToby Isaac   } else {
439566063dSJacob Faibussowitsch     PetscCall(PetscParallelSortInt_Bitonic_Merge(comm, tag, mid, rankEnd, rank, n, keys, buffer, forward));
44ce605777SToby Isaac   }
45ce605777SToby Isaac   PetscFunctionReturn(0);
46ce605777SToby Isaac }
47ce605777SToby Isaac 
48ce605777SToby Isaac /* This is the bitonic sort that works on non-power-of-2 sizes found at http://www.iti.fh-flensburg.de/lang/algorithmen/sortieren/bitonic/oddn.htm */
499371c9d4SSatish Balay static PetscErrorCode PetscParallelSortInt_Bitonic_Recursive(MPI_Comm comm, PetscMPIInt tag, PetscMPIInt rankStart, PetscMPIInt rankEnd, PetscMPIInt rank, PetscMPIInt n, PetscInt keys[], PetscInt buffer[], PetscBool forward) {
50ce605777SToby Isaac   PetscInt diff;
51ce605777SToby Isaac   PetscInt mid;
52ce605777SToby Isaac 
53ce605777SToby Isaac   PetscFunctionBegin;
54ce605777SToby Isaac   diff = rankEnd - rankStart;
55ce605777SToby Isaac   if (diff <= 0) PetscFunctionReturn(0);
56ce605777SToby Isaac   if (diff == 1) {
57ce605777SToby Isaac     if (forward) {
589566063dSJacob Faibussowitsch       PetscCall(PetscSortInt(n, keys));
59ce605777SToby Isaac     } else {
609566063dSJacob Faibussowitsch       PetscCall(PetscSortReverseInt(n, keys));
61ce605777SToby Isaac     }
62ce605777SToby Isaac     PetscFunctionReturn(0);
63ce605777SToby Isaac   }
64ce605777SToby Isaac   mid = rankStart + diff / 2;
65ce605777SToby Isaac   /* divide and conquer */
66ce605777SToby Isaac   if (rank < mid) {
679566063dSJacob Faibussowitsch     PetscCall(PetscParallelSortInt_Bitonic_Recursive(comm, tag, rankStart, mid, rank, n, keys, buffer, (PetscBool)!forward));
68ce605777SToby Isaac   } else {
699566063dSJacob Faibussowitsch     PetscCall(PetscParallelSortInt_Bitonic_Recursive(comm, tag, mid, rankEnd, rank, n, keys, buffer, forward));
70ce605777SToby Isaac   }
71ce605777SToby Isaac   /* bitonic merge */
729566063dSJacob Faibussowitsch   PetscCall(PetscParallelSortInt_Bitonic_Merge(comm, tag, rankStart, rankEnd, rank, n, keys, buffer, forward));
73ce605777SToby Isaac   PetscFunctionReturn(0);
74ce605777SToby Isaac }
75ce605777SToby Isaac 
769371c9d4SSatish Balay static PetscErrorCode PetscParallelSortInt_Bitonic(MPI_Comm comm, PetscInt n, PetscInt keys[]) {
77ce605777SToby Isaac   PetscMPIInt size, rank, tag, mpin;
78ce605777SToby Isaac   PetscInt   *buffer;
79ce605777SToby Isaac 
80ce605777SToby Isaac   PetscFunctionBegin;
81ce605777SToby Isaac   PetscValidIntPointer(keys, 3);
829566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(comm, &tag));
839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
859566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(n, &mpin));
869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &buffer));
879566063dSJacob Faibussowitsch   PetscCall(PetscParallelSortInt_Bitonic_Recursive(comm, tag, 0, size, rank, mpin, keys, buffer, PETSC_TRUE));
889566063dSJacob Faibussowitsch   PetscCall(PetscFree(buffer));
89ce605777SToby Isaac   PetscFunctionReturn(0);
90ce605777SToby Isaac }
91ce605777SToby Isaac 
929371c9d4SSatish Balay static PetscErrorCode PetscParallelSampleSelect(PetscLayout mapin, PetscLayout mapout, PetscInt keysin[], PetscInt *outpivots[]) {
93ce605777SToby Isaac   PetscMPIInt  size, rank;
94ce605777SToby Isaac   PetscInt    *pivots, *finalpivots, i;
95ce605777SToby Isaac   PetscInt     non_empty, my_first, count;
96ce605777SToby Isaac   PetscMPIInt *keys_per, max_keys_per;
97ce605777SToby Isaac 
98ce605777SToby Isaac   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(mapin->comm, &size));
1009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(mapin->comm, &rank));
101ce605777SToby Isaac 
102ce605777SToby Isaac   /* Choose P - 1 pivots that would be ideal for the distribution on this process */
1039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size - 1, &pivots));
104ce605777SToby Isaac   for (i = 0; i < size - 1; i++) {
105ce605777SToby Isaac     PetscInt index;
106ce605777SToby Isaac 
107ce605777SToby Isaac     if (!mapin->n) {
108ce605777SToby Isaac       /* if this rank is empty, put "infinity" in the list.  each process knows how many empty
109ce605777SToby Isaac        * processes are in the layout, so those values will be ignored from the end of the sorted
110ce605777SToby Isaac        * pivots */
111ce605777SToby Isaac       pivots[i] = PETSC_MAX_INT;
112ce605777SToby Isaac     } else {
113ce605777SToby Isaac       /* match the distribution to the desired output described by mapout by getting the index
114ce605777SToby Isaac        * that is approximately the appropriate fraction through the list */
115ce605777SToby Isaac       index     = ((PetscReal)mapout->range[i + 1]) * ((PetscReal)mapin->n) / ((PetscReal)mapout->N);
116ce605777SToby Isaac       index     = PetscMin(index, (mapin->n - 1));
117ce605777SToby Isaac       index     = PetscMax(index, 0);
118ce605777SToby Isaac       pivots[i] = keysin[index];
119ce605777SToby Isaac     }
120ce605777SToby Isaac   }
121ce605777SToby Isaac   /* sort the pivots in parallel */
1229566063dSJacob Faibussowitsch   PetscCall(PetscParallelSortInt_Bitonic(mapin->comm, size - 1, pivots));
12376bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
124ce605777SToby Isaac     PetscBool sorted;
125ce605777SToby Isaac 
1269566063dSJacob Faibussowitsch     PetscCall(PetscParallelSortedInt(mapin->comm, size - 1, pivots, &sorted));
12728b400f6SJacob Faibussowitsch     PetscCheck(sorted, mapin->comm, PETSC_ERR_PLIB, "bitonic sort failed");
128ce605777SToby Isaac   }
129ce605777SToby Isaac 
130ce605777SToby Isaac   /* if there are Z nonempty processes, we have (P - 1) * Z real pivots, and we want to select
131ce605777SToby Isaac    * at indices Z - 1, 2*Z - 1, ... (P - 1) * Z - 1 */
132ce605777SToby Isaac   non_empty = size;
1339371c9d4SSatish Balay   for (i = 0; i < size; i++)
1349371c9d4SSatish Balay     if (mapout->range[i] == mapout->range[i + 1]) non_empty--;
1359566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(size + 1, &keys_per));
136ce605777SToby Isaac   my_first = -1;
137ce605777SToby Isaac   if (non_empty) {
138ce605777SToby Isaac     for (i = 0; i < size - 1; i++) {
139ce605777SToby Isaac       size_t sample      = (i + 1) * non_empty - 1;
140ce605777SToby Isaac       size_t sample_rank = sample / (size - 1);
141ce605777SToby Isaac 
142ce605777SToby Isaac       keys_per[sample_rank]++;
1439371c9d4SSatish Balay       if (my_first < 0 && (PetscMPIInt)sample_rank == rank) { my_first = (PetscInt)(sample - sample_rank * (size - 1)); }
144ce605777SToby Isaac     }
145ce605777SToby Isaac   }
146ce605777SToby Isaac   for (i = 0, max_keys_per = 0; i < size; i++) max_keys_per = PetscMax(keys_per[i], max_keys_per);
1479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size * max_keys_per, &finalpivots));
148ce605777SToby Isaac   /* now that we know how many pivots each process will provide, gather the selected pivots at the start of the array
149ce605777SToby Isaac    * and allgather them */
150ce605777SToby Isaac   for (i = 0; i < keys_per[rank]; i++) pivots[i] = pivots[my_first + i * non_empty];
151ce605777SToby Isaac   for (i = keys_per[rank]; i < max_keys_per; i++) pivots[i] = PETSC_MAX_INT;
1529566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(pivots, max_keys_per, MPIU_INT, finalpivots, max_keys_per, MPIU_INT, mapin->comm));
153ce605777SToby Isaac   for (i = 0, count = 0; i < size; i++) {
154ce605777SToby Isaac     PetscInt j;
155ce605777SToby Isaac 
156ce605777SToby Isaac     for (j = 0; j < max_keys_per; j++) {
1579371c9d4SSatish Balay       if (j < keys_per[i]) { finalpivots[count++] = finalpivots[i * max_keys_per + j]; }
158ce605777SToby Isaac     }
159ce605777SToby Isaac   }
160ce605777SToby Isaac   *outpivots = finalpivots;
1619566063dSJacob Faibussowitsch   PetscCall(PetscFree(keys_per));
1629566063dSJacob Faibussowitsch   PetscCall(PetscFree(pivots));
163ce605777SToby Isaac   PetscFunctionReturn(0);
164ce605777SToby Isaac }
165ce605777SToby Isaac 
1669371c9d4SSatish Balay static PetscErrorCode PetscParallelRedistribute(PetscLayout map, PetscInt n, PetscInt arrayin[], PetscInt arrayout[]) {
167ce605777SToby Isaac   PetscMPIInt  size, rank;
168ce605777SToby Isaac   PetscInt     myOffset, nextOffset;
169ce605777SToby Isaac   PetscInt     i;
170ce605777SToby Isaac   PetscMPIInt  total, filled;
171ce605777SToby Isaac   PetscMPIInt  sender, nfirst, nsecond;
172ce605777SToby Isaac   PetscMPIInt  firsttag, secondtag;
173ce605777SToby Isaac   MPI_Request  firstreqrcv;
174ce605777SToby Isaac   MPI_Request *firstreqs;
175ce605777SToby Isaac   MPI_Request *secondreqs;
176ce605777SToby Isaac   MPI_Status   firststatus;
177ce605777SToby Isaac 
178ce605777SToby Isaac   PetscFunctionBegin;
1799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(map->comm, &size));
1809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(map->comm, &rank));
1819566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(map->comm, &firsttag));
1829566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(map->comm, &secondtag));
183ce605777SToby Isaac   myOffset = 0;
1849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size, &firstreqs, size, &secondreqs));
1859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Scan(&n, &nextOffset, 1, MPIU_INT, MPI_SUM, map->comm));
186ce605777SToby Isaac   myOffset = nextOffset - n;
187ce605777SToby Isaac   total    = map->range[rank + 1] - map->range[rank];
188*48a46eb9SPierre Jolivet   if (total > 0) PetscCallMPI(MPI_Irecv(arrayout, total, MPIU_INT, MPI_ANY_SOURCE, firsttag, map->comm, &firstreqrcv));
189ce605777SToby Isaac   for (i = 0, nsecond = 0, nfirst = 0; i < size; i++) {
190ce605777SToby Isaac     PetscInt itotal;
191ce605777SToby Isaac     PetscInt overlap, oStart, oEnd;
192ce605777SToby Isaac 
193ce605777SToby Isaac     itotal = map->range[i + 1] - map->range[i];
194ce605777SToby Isaac     if (itotal <= 0) continue;
195ce605777SToby Isaac     oStart  = PetscMax(myOffset, map->range[i]);
196ce605777SToby Isaac     oEnd    = PetscMin(nextOffset, map->range[i + 1]);
197ce605777SToby Isaac     overlap = oEnd - oStart;
198ce605777SToby Isaac     if (map->range[i] >= myOffset && map->range[i] < nextOffset) {
199ce605777SToby Isaac       /* send first message */
2009566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(&arrayin[map->range[i] - myOffset], overlap, MPIU_INT, i, firsttag, map->comm, &(firstreqs[nfirst++])));
201ce605777SToby Isaac     } else if (overlap > 0) {
202ce605777SToby Isaac       /* send second message */
2039566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(&arrayin[oStart - myOffset], overlap, MPIU_INT, i, secondtag, map->comm, &(secondreqs[nsecond++])));
204ce605777SToby Isaac     } else if (overlap == 0 && myOffset > map->range[i] && myOffset < map->range[i + 1]) {
205ce605777SToby Isaac       /* send empty second message */
2069566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(&arrayin[oStart - myOffset], 0, MPIU_INT, i, secondtag, map->comm, &(secondreqs[nsecond++])));
207ce605777SToby Isaac     }
208ce605777SToby Isaac   }
209ce605777SToby Isaac   filled = 0;
210ce605777SToby Isaac   sender = -1;
211ce605777SToby Isaac   if (total > 0) {
2129566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Wait(&firstreqrcv, &firststatus));
213ce605777SToby Isaac     sender = firststatus.MPI_SOURCE;
2149566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_count(&firststatus, MPIU_INT, &filled));
215ce605777SToby Isaac   }
216ce605777SToby Isaac   while (filled < total) {
217ce605777SToby Isaac     PetscMPIInt mfilled;
218ce605777SToby Isaac     MPI_Status  stat;
219ce605777SToby Isaac 
220ce605777SToby Isaac     sender++;
2219566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Recv(&arrayout[filled], total - filled, MPIU_INT, sender, secondtag, map->comm, &stat));
2229566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_count(&stat, MPIU_INT, &mfilled));
223ce605777SToby Isaac     filled += mfilled;
224ce605777SToby Isaac   }
2259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(nfirst, firstreqs, MPI_STATUSES_IGNORE));
2269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(nsecond, secondreqs, MPI_STATUSES_IGNORE));
2279566063dSJacob Faibussowitsch   PetscCall(PetscFree2(firstreqs, secondreqs));
228ce605777SToby Isaac   PetscFunctionReturn(0);
229ce605777SToby Isaac }
230ce605777SToby Isaac 
2319371c9d4SSatish Balay static PetscErrorCode PetscParallelSortInt_Samplesort(PetscLayout mapin, PetscLayout mapout, PetscInt keysin[], PetscInt keysout[]) {
232ce605777SToby Isaac   PetscMPIInt  size, rank;
2336df2664eSToby Isaac   PetscInt    *pivots = NULL, *buffer;
234ce605777SToby Isaac   PetscInt     i, j;
235ce605777SToby Isaac   PetscMPIInt *keys_per_snd, *keys_per_rcv, *offsets_snd, *offsets_rcv, nrecv;
236ce605777SToby Isaac 
237ce605777SToby Isaac   PetscFunctionBegin;
2389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(mapin->comm, &size));
2399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(mapin->comm, &rank));
2409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(size, &keys_per_snd, size, &keys_per_rcv, size + 1, &offsets_snd, size + 1, &offsets_rcv));
241ce605777SToby Isaac   /* sort locally */
2429566063dSJacob Faibussowitsch   PetscCall(PetscSortInt(mapin->n, keysin));
243ce605777SToby Isaac   /* get P - 1 pivots */
2449566063dSJacob Faibussowitsch   PetscCall(PetscParallelSampleSelect(mapin, mapout, keysin, &pivots));
245ce605777SToby Isaac   /* determine which entries in the sorted array go to which other processes based on the pivots */
246ce605777SToby Isaac   for (i = 0, j = 0; i < size - 1; i++) {
247ce605777SToby Isaac     PetscInt prev = j;
248ce605777SToby Isaac 
249ce605777SToby Isaac     while ((j < mapin->n) && (keysin[j] < pivots[i])) j++;
250ce605777SToby Isaac     offsets_snd[i]  = prev;
251ce605777SToby Isaac     keys_per_snd[i] = j - prev;
252ce605777SToby Isaac   }
253ce605777SToby Isaac   offsets_snd[size - 1]  = j;
254ce605777SToby Isaac   keys_per_snd[size - 1] = mapin->n - j;
255ce605777SToby Isaac   offsets_snd[size]      = mapin->n;
256ce605777SToby Isaac   /* get the incoming sizes */
2579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Alltoall(keys_per_snd, 1, MPI_INT, keys_per_rcv, 1, MPI_INT, mapin->comm));
258ce605777SToby Isaac   offsets_rcv[0] = 0;
2599371c9d4SSatish Balay   for (i = 0; i < size; i++) { offsets_rcv[i + 1] = offsets_rcv[i] + keys_per_rcv[i]; }
260ce605777SToby Isaac   nrecv = offsets_rcv[size];
261ce605777SToby Isaac   /* all to all exchange */
2629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecv, &buffer));
2639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Alltoallv(keysin, keys_per_snd, offsets_snd, MPIU_INT, buffer, keys_per_rcv, offsets_rcv, MPIU_INT, mapin->comm));
2649566063dSJacob Faibussowitsch   PetscCall(PetscFree(pivots));
2659566063dSJacob Faibussowitsch   PetscCall(PetscFree4(keys_per_snd, keys_per_rcv, offsets_snd, offsets_rcv));
266ce605777SToby Isaac 
267ce605777SToby Isaac   /* local sort */
2689566063dSJacob Faibussowitsch   PetscCall(PetscSortInt(nrecv, buffer));
269ce605777SToby Isaac #if defined(PETSC_USE_DEBUG)
270ce605777SToby Isaac   {
271ce605777SToby Isaac     PetscBool sorted;
272ce605777SToby Isaac 
2739566063dSJacob Faibussowitsch     PetscCall(PetscParallelSortedInt(mapin->comm, nrecv, buffer, &sorted));
27428b400f6SJacob Faibussowitsch     PetscCheck(sorted, mapin->comm, PETSC_ERR_PLIB, "samplesort (pre-redistribute) sort failed");
275ce605777SToby Isaac   }
276ce605777SToby Isaac #endif
277ce605777SToby Isaac 
278ce605777SToby Isaac   /* redistribute to the desired order */
2799566063dSJacob Faibussowitsch   PetscCall(PetscParallelRedistribute(mapout, nrecv, buffer, keysout));
2809566063dSJacob Faibussowitsch   PetscCall(PetscFree(buffer));
281ce605777SToby Isaac   PetscFunctionReturn(0);
282ce605777SToby Isaac }
283ce605777SToby Isaac 
284ce605777SToby Isaac /*@
285064ec1f3Sprj-   PetscParallelSortInt - Globally sort a distributed array of integers
286ce605777SToby Isaac 
287ce605777SToby Isaac   Collective
288ce605777SToby Isaac 
289ce605777SToby Isaac   Input Parameters:
290ce605777SToby Isaac + mapin - PetscLayout describing the distribution of the input keys
291ce605777SToby Isaac . mapout - PetscLayout describing the distribution of the output keys
292ce605777SToby Isaac - keysin - the pre-sorted array of integers
293ce605777SToby Isaac 
294064ec1f3Sprj-   Output Parameter:
295ce605777SToby Isaac . keysout - the array in which the sorted integers will be stored.  If mapin == mapout, then keysin may be equal to keysout.
296ce605777SToby Isaac 
297ce605777SToby Isaac   Level: developer
298ce605777SToby Isaac 
299ce605777SToby Isaac   Notes: This implements a distributed samplesort, which, with local array sizes n_in and n_out, global size N, and global number of processes P, does:
300ce605777SToby Isaac 
301ce605777SToby Isaac   - sorts locally
302ce605777SToby Isaac   - chooses pivots by sorting (in parallel) (P-1) pivot suggestions from each process using bitonic sort and allgathering a subset of (P-1) of those
303ce605777SToby Isaac   - using to the pivots to repartition the keys by all-to-all exchange
304ce605777SToby Isaac   - sorting the repartitioned keys locally (the array is now globally sorted, but does not match the mapout layout)
305ce605777SToby Isaac   - redistributing to match the mapout layout
306ce605777SToby Isaac 
307ce605777SToby Isaac   If keysin != keysout, then keysin will not be changed during PetscParallelSortInt.
308ce605777SToby Isaac 
309db781477SPatrick Sanan .seealso: `PetscParallelSortedInt()`
310ce605777SToby Isaac @*/
3119371c9d4SSatish Balay PetscErrorCode PetscParallelSortInt(PetscLayout mapin, PetscLayout mapout, PetscInt keysin[], PetscInt keysout[]) {
312ce605777SToby Isaac   PetscMPIInt size;
313ce605777SToby Isaac   PetscMPIInt result;
314ce605777SToby Isaac   PetscInt   *keysincopy = NULL;
315ce605777SToby Isaac 
316ce605777SToby Isaac   PetscFunctionBegin;
317ce605777SToby Isaac   PetscValidPointer(mapin, 1);
318ce605777SToby Isaac   PetscValidPointer(mapout, 2);
3199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_compare(mapin->comm, mapout->comm, &result));
320c9cc58a2SBarry Smith   PetscCheck(result == MPI_IDENT || result == MPI_CONGRUENT, mapin->comm, PETSC_ERR_ARG_NOTSAMECOMM, "layouts are not on the same communicator");
3219566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(mapin));
3229566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(mapout));
323ce605777SToby Isaac   if (mapin->n) PetscValidIntPointer(keysin, 3);
324ce605777SToby Isaac   if (mapout->n) PetscValidIntPointer(keysout, 4);
32508401ef6SPierre Jolivet   PetscCheck(mapin->N == mapout->N, mapin->comm, PETSC_ERR_ARG_SIZ, "Input and output layouts have different global sizes (%" PetscInt_FMT " != %" PetscInt_FMT ")", mapin->N, mapout->N);
3269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(mapin->comm, &size));
327ce605777SToby Isaac   if (size == 1) {
328*48a46eb9SPierre Jolivet     if (keysout != keysin) PetscCall(PetscMemcpy(keysout, keysin, mapin->n * sizeof(PetscInt)));
3299566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(mapout->n, keysout));
330ce605777SToby Isaac     if (size == 1) PetscFunctionReturn(0);
331ce605777SToby Isaac   }
332ce605777SToby Isaac   if (keysout != keysin) {
3339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mapin->n, &keysincopy));
3349566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(keysincopy, keysin, mapin->n * sizeof(PetscInt)));
335ce605777SToby Isaac     keysin = keysincopy;
336ce605777SToby Isaac   }
3379566063dSJacob Faibussowitsch   PetscCall(PetscParallelSortInt_Samplesort(mapin, mapout, keysin, keysout));
338ce605777SToby Isaac #if defined(PETSC_USE_DEBUG)
339ce605777SToby Isaac   {
340ce605777SToby Isaac     PetscBool sorted;
341ce605777SToby Isaac 
3429566063dSJacob Faibussowitsch     PetscCall(PetscParallelSortedInt(mapout->comm, mapout->n, keysout, &sorted));
34328b400f6SJacob Faibussowitsch     PetscCheck(sorted, mapout->comm, PETSC_ERR_PLIB, "samplesort sort failed");
344ce605777SToby Isaac   }
345ce605777SToby Isaac #endif
3469566063dSJacob Faibussowitsch   PetscCall(PetscFree(keysincopy));
347ce605777SToby Isaac   PetscFunctionReturn(0);
348ce605777SToby Isaac }
349