MP-Gadget  5.0.1.dev1-76bc7d4726-dirty
Classes | Macros | Typedefs | Functions | Variables
mpsort.c File Reference
#include <stdio.h>
#include <stddef.h>
#include <stdint.h>
#include <stdlib.h>
#include <time.h>
#include <string.h>
#include <mpi.h>
#include "mpsort.h"
#include "system.h"
#include "mymalloc.h"
#include "openmpsort.h"
#include "endrun.h"
Include dependency graph for mpsort.c:

Go to the source code of this file.

Classes

struct  crstruct
 
struct  piter
 
struct  crmpistruct
 
struct  TIMER
 
struct  TIMERS
 
struct  SegmentGroupDescr
 

Macros

#define DEFTYPE(type)
 

Typedefs

typedef int(* _compar_fn_t) (const void *r1, const void *r2, size_t rsize)
 
typedef void(* _bisect_fn_t) (void *r, const void *r1, const void *r2, size_t rsize)
 

Functions

static int _compar_radix (const void *r1, const void *r2, size_t rsize, int dir)
 
static int _compar_radix_u8 (const void *r1, const void *r2, size_t rsize, int dir)
 
static int _compar_radix_le (const void *r1, const void *r2, size_t rsize)
 
static int _compar_radix_be (const void *r1, const void *r2, size_t rsize)
 
static int _compar_radix_le_u8 (const void *r1, const void *r2, size_t rsize)
 
static int _compar_radix_be_u8 (const void *r1, const void *r2, size_t rsize)
 
static void _bisect_radix (void *r, const void *r1, const void *r2, size_t rsize, int dir)
 
static void _bisect_radix_le (void *r, const void *r1, const void *r2, size_t rsize)
 
static void _bisect_radix_be (void *r, const void *r1, const void *r2, size_t rsize)
 
void _setup_radix_sort (struct crstruct *d, void *base, size_t nmemb, size_t size, void(*radix)(const void *ptr, void *radix, void *arg), size_t rsize, void *arg)
 
static int _compute_and_compar_radix (const void *p1, const void *p2)
 
static void radix_sort (void *base, size_t nmemb, size_t size, void(*radix)(const void *ptr, void *radix, void *arg), size_t rsize, void *arg)
 
static ptrdiff_t _bsearch_last_lt (void *P, void *base, size_t nmemb, struct crstruct *d)
 
static ptrdiff_t _bsearch_last_le (void *P, void *base, size_t nmemb, struct crstruct *d)
 
static void _histogram (char *P, int Plength, void *mybase, size_t mynmemb, ptrdiff_t *myCLT, ptrdiff_t *myCLE, struct crstruct *d)
 
static void piter_init (struct piter *pi, char *Pmin, char *Pmax, int Plength, struct crstruct *d)
 
static void piter_destroy (struct piter *pi)
 
static void piter_bisect (struct piter *pi, char *P)
 
static int piter_all_done (struct piter *pi)
 
static void piter_accept (struct piter *pi, char *P, ptrdiff_t *C, ptrdiff_t *CLT, ptrdiff_t *CLE)
 
static void _setup_mpsort_mpi (struct crmpistruct *o, struct crstruct *d, void *myoutbase, size_t myoutnmemb, MPI_Comm comm)
 
static void _destroy_mpsort_mpi (struct crmpistruct *o)
 
static void _find_Pmax_Pmin_C (void *mybase, size_t mynmemb, size_t myoutnmemb, char *Pmax, char *Pmin, ptrdiff_t *C, struct crstruct *d, struct crmpistruct *o)
 
static int _solve_for_layout_mpi (int NTask, ptrdiff_t *C, ptrdiff_t *myT_CLT, ptrdiff_t *myT_CLE, ptrdiff_t *myT_C, MPI_Comm comm)
 
static int _assign_colors (size_t glocalsize, size_t *sizes, int *ncolor, MPI_Comm comm)
 
static size_t _collect_sizes (size_t localsize, size_t *sizes, size_t *myoffset, MPI_Comm comm)
 
static int MPIU_GetLoc (const void *base, MPI_Datatype type, MPI_Op op, MPI_Comm comm)
 
static void _create_segment_group (struct SegmentGroupDescr *descr, size_t *sizes, size_t avgsegsize, int Ngroup, MPI_Comm comm)
 
static void _destroy_segment_group (struct SegmentGroupDescr *descr)
 
static void mpsort_increment_timer (const char *name, int erase)
 
void mpsort_setup_timers (int ntimers)
 
void mpsort_free_timers (void)
 
void mpsort_mpi_report_last_run ()
 
int mpsort_mpi_find_ntimers (struct TIMERS *timers)
 
void mpsort_mpi_impl (void *mybase, size_t mynmemb, size_t size, void(*radix)(const void *ptr, void *radix, void *arg), size_t rsize, void *arg, MPI_Comm comm, const int line, const char *file)
 
static int mpsort_mpi_histogram_sort (struct crstruct d, struct crmpistruct o)
 
static void MPIU_Scatter (MPI_Comm comm, int root, const void *sendbuffer, void *recvbuffer, int nrecv, size_t elsize, int *totalnsend)
 
static void MPIU_Gather (MPI_Comm comm, int root, const void *sendbuffer, void *recvbuffer, int nsend, size_t elsize, int *totalnrecv)
 
static uint64_t checksum (void *base, size_t nbytes, MPI_Comm comm)
 
void mpsort_mpi_newarray_impl (void *mybase, size_t mynmemb, void *myoutbase, size_t myoutnmemb, size_t elsize, void(*radix)(const void *ptr, void *radix, void *arg), size_t rsize, void *arg, MPI_Comm comm, const int line, const char *file)
 
static void _mpsort_mpi_parse_env ()
 
void mpsort_mpi_set_options (int options)
 
int mpsort_mpi_has_options (int options)
 
void mpsort_mpi_unset_options (int options)
 

Variables

static struct crstruct _cacr_d
 
static int _mpsort_mpi_options = 0
 
static MPI_Datatype MPI_TYPE_PTRDIFF = 0
 
static struct TIMERS _TIMERS
 

Macro Definition Documentation

◆ DEFTYPE

#define DEFTYPE (   type)
Value:
static int _compar_radix_ ## type ( \
const type * u1, \
const type * u2, \
size_t junk) { \
return (signed) (*u1 > *u2) - (signed) (*u1 < *u2); \
} \
static void _bisect_radix_ ## type ( \
type * u, \
const type * u1, \
const type * u2, \
size_t junk) { \
*u = *u1 + ((*u2 - *u1) >> 1); \
}

Definition at line 31 of file mpsort.c.

Typedef Documentation

◆ _bisect_fn_t

typedef void(* _bisect_fn_t) (void *r, const void *r1, const void *r2, size_t rsize)

Definition at line 18 of file mpsort.c.

◆ _compar_fn_t

typedef int(* _compar_fn_t) (const void *r1, const void *r2, size_t rsize)

Definition at line 17 of file mpsort.c.

Function Documentation

◆ _assign_colors()

static int _assign_colors ( size_t  glocalsize,
size_t *  sizes,
int *  ncolor,
MPI_Comm  comm 
)
static

Definition at line 555 of file mpsort.c.

556 {
557  int NTask;
558  int ThisTask;
559  MPI_Comm_rank(comm, &ThisTask);
560  MPI_Comm_size(comm, &NTask);
561 
562  int i;
563  int mycolor = -1;
564  size_t current_size = 0;
565  size_t current_outsize = 0;
566  int current_color = 0;
567  int lastcolor = 0;
568  for(i = 0; i < NTask; i ++) {
569  current_size += sizes[i];
570 
571  lastcolor = current_color;
572 
573  if(i == ThisTask) {
574  mycolor = lastcolor;
575  }
576 
577  if(current_size > glocalsize || current_outsize > glocalsize) {
578  current_size = 0;
579  current_outsize = 0;
580  current_color ++;
581  }
582  }
583  /* no data for color of -1; exclude them later with special cases */
584  if(sizes[ThisTask] == 0) {
585  mycolor = -1;
586  }
587 
588  *ncolor = lastcolor + 1;
589  return mycolor;
590 }
int ThisTask
Definition: test_exchange.c:23
int NTask
Definition: test_exchange.c:23

References NTask, and ThisTask.

Referenced by _create_segment_group().

Here is the caller graph for this function:

◆ _bisect_radix()

static void _bisect_radix ( void *  r,
const void *  r1,
const void *  r2,
size_t  rsize,
int  dir 
)
static

Definition at line 94 of file mpsort.c.

94  {
95  size_t i;
96  const unsigned char * u1 = (const unsigned char *) r1;
97  const unsigned char * u2 = (const unsigned char *) r2;
98  unsigned char * u = (unsigned char *) r;
99  unsigned int carry = 0;
100  if(dir > 0) {
101  u1 += rsize - 1;
102  u2 += rsize - 1;
103  }
104  /* from most significant */
105  for(i = 0; i < rsize; i ++) {
106  unsigned int tmp = (unsigned int) *u2 + *u1 + carry;
107  if(tmp >= 256) carry = 1;
108  else carry = 0;
109  *u = tmp % (UINT8_MAX+1);
110  u -= dir;
111  u1 -= dir;
112  u2 -= dir;
113  }
114  u += dir;
115  for(i = 0; i < rsize; i ++) {
116  unsigned int tmp = *u + carry * 256;
117  carry = tmp & 1;
118  *u = (tmp >> 1) ;
119  u += dir;
120  }
121 }

Referenced by _bisect_radix_be(), and _bisect_radix_le().

Here is the caller graph for this function:

◆ _bisect_radix_be()

static void _bisect_radix_be ( void *  r,
const void *  r1,
const void *  r2,
size_t  rsize 
)
static

Definition at line 125 of file mpsort.c.

125  {
126  _bisect_radix(r, r1, r2, rsize, +1);
127 }
static void _bisect_radix(void *r, const void *r1, const void *r2, size_t rsize, int dir)
Definition: mpsort.c:94

References _bisect_radix().

Here is the call graph for this function:

◆ _bisect_radix_le()

static void _bisect_radix_le ( void *  r,
const void *  r1,
const void *  r2,
size_t  rsize 
)
static

Definition at line 122 of file mpsort.c.

122  {
123  _bisect_radix(r, r1, r2, rsize, -1);
124 }

References _bisect_radix().

Here is the call graph for this function:

◆ _bsearch_last_le()

static ptrdiff_t _bsearch_last_le ( void *  P,
void *  base,
size_t  nmemb,
struct crstruct d 
)
static

Definition at line 257 of file mpsort.c.

259  {
260 
261  if (nmemb == 0) return -1;
262 
263  char tmpradix[d->rsize];
264  ptrdiff_t left = 0;
265  ptrdiff_t right = nmemb - 1;
266 
267  d->radix((char*) base, tmpradix, d->arg);
268  if(d->compar(tmpradix, P, d->rsize) > 0) {
269  return -1;
270  }
271  d->radix((char*) base + right * d->size, tmpradix, d->arg);
272  if(d->compar(tmpradix, P, d->rsize) <= 0) {
273  return nmemb - 1;
274  }
275 
276  /* left <= i <= right*/
277  /* [left] <= P < [right] */
278  while(right > left + 1) {
279  ptrdiff_t mid = ((right - left + 1) >> 1) + left;
280  d->radix((char*) base + mid * d->size, tmpradix, d->arg);
281  /* if [mid] <= P , move left to mid */
282  /* if [mid] > P , move right to mid*/
283  int c1 = d->compar(tmpradix, P, d->rsize);
284  if(c1 <= 0) {
285  left = mid;
286  } else {
287  right = mid;
288  }
289  }
290  return left;
291 }
#define P
Definition: partmanager.h:88
size_t size
Definition: mpsort.c:23
_compar_fn_t compar
Definition: mpsort.c:27
void(* radix)(const void *ptr, void *radix, void *arg)
Definition: mpsort.c:26
void * arg
Definition: mpsort.c:25
size_t rsize
Definition: mpsort.c:24

References crstruct::arg, crstruct::base, crstruct::compar, crstruct::nmemb, P, crstruct::radix, crstruct::rsize, and crstruct::size.

Referenced by _histogram().

Here is the caller graph for this function:

◆ _bsearch_last_lt()

static ptrdiff_t _bsearch_last_lt ( void *  P,
void *  base,
size_t  nmemb,
struct crstruct d 
)
static

Definition at line 216 of file mpsort.c.

218  {
219 
220  if (nmemb == 0) return -1;
221 
222  char tmpradix[d->rsize];
223  ptrdiff_t left = 0;
224  ptrdiff_t right = nmemb - 1;
225 
226  d->radix((char*) base, tmpradix, d->arg);
227  if(d->compar(tmpradix, P, d->rsize) >= 0) {
228  return - 1;
229  }
230  d->radix((char*) base + right * d->size, tmpradix, d->arg);
231  if(d->compar(tmpradix, P, d->rsize) < 0) {
232  return nmemb - 1;
233  }
234 
235  /* left <= i <= right*/
236  /* [left] < P <= [right] */
237  while(right > left + 1) {
238  ptrdiff_t mid = ((right - left + 1) >> 1) + left;
239  d->radix((char*) base + mid * d->size, tmpradix, d->arg);
240  /* if [mid] < P , move left to mid */
241  /* if [mid] >= P , move right to mid */
242  int c1 = d->compar(tmpradix, P, d->rsize);
243  if(c1 < 0) {
244  left = mid;
245  } else {
246  right = mid;
247  }
248  }
249  return left;
250 }

References crstruct::arg, crstruct::base, crstruct::compar, crstruct::nmemb, P, crstruct::radix, crstruct::rsize, and crstruct::size.

Referenced by _histogram().

Here is the caller graph for this function:

◆ _collect_sizes()

static size_t _collect_sizes ( size_t  localsize,
size_t *  sizes,
size_t *  myoffset,
MPI_Comm  comm 
)
static

Definition at line 593 of file mpsort.c.

594 {
595 
596  int ThisTask, NTask;
597 
598  MPI_Comm_size(comm, &NTask);
599  MPI_Comm_rank(comm, &ThisTask);
600 
601  size_t totalsize;
602 
603  sizes[ThisTask] = localsize;
604 
605  MPI_Allreduce(&sizes[ThisTask], &totalsize, 1, MPI_TYPE_PTRDIFF, MPI_SUM, comm);
606  MPI_Allgather(MPI_IN_PLACE, 0, MPI_DATATYPE_NULL, sizes, 1, MPI_TYPE_PTRDIFF, comm);
607 
608  int i;
609  *myoffset = 0;
610  for(i = 0; i < ThisTask; i ++) {
611  (*myoffset) += sizes[i];
612  }
613 
614  return totalsize;
615 }
static MPI_Datatype MPI_TYPE_PTRDIFF
Definition: mpsort.c:476

References MPI_TYPE_PTRDIFF, NTask, and ThisTask.

Referenced by mpsort_mpi_newarray_impl().

Here is the caller graph for this function:

◆ _compar_radix()

static int _compar_radix ( const void *  r1,
const void *  r2,
size_t  rsize,
int  dir 
)
static

Definition at line 48 of file mpsort.c.

48  {
49  size_t i;
50  /* from most significant */
51  const unsigned char * u1 = (const unsigned char *) r1;
52  const unsigned char * u2 = (const unsigned char *) r2;
53  if(dir < 0) {
54  u1 += rsize - 1;
55  u2 += rsize - 1;;
56  }
57  for(i = 0; i < rsize; i ++) {
58  if(*u1 < *u2) return -1;
59  if(*u1 > *u2) return 1;
60  u1 += dir;
61  u2 += dir;
62  }
63  return 0;
64 }

Referenced by _compar_radix_be(), and _compar_radix_le().

Here is the caller graph for this function:

◆ _compar_radix_be()

static int _compar_radix_be ( const void *  r1,
const void *  r2,
size_t  rsize 
)
static

Definition at line 85 of file mpsort.c.

85  {
86  return _compar_radix(r1, r2, rsize, +1);
87 }
static int _compar_radix(const void *r1, const void *r2, size_t rsize, int dir)
Definition: mpsort.c:48

References _compar_radix().

Here is the call graph for this function:

◆ _compar_radix_be_u8()

static int _compar_radix_be_u8 ( const void *  r1,
const void *  r2,
size_t  rsize 
)
static

Definition at line 91 of file mpsort.c.

91  {
92  return _compar_radix_u8(r1, r2, rsize, +1);
93 }
static int _compar_radix_u8(const void *r1, const void *r2, size_t rsize, int dir)
Definition: mpsort.c:65

References _compar_radix_u8().

Here is the call graph for this function:

◆ _compar_radix_le()

static int _compar_radix_le ( const void *  r1,
const void *  r2,
size_t  rsize 
)
static

Definition at line 82 of file mpsort.c.

82  {
83  return _compar_radix(r1, r2, rsize, -1);
84 }

References _compar_radix().

Here is the call graph for this function:

◆ _compar_radix_le_u8()

static int _compar_radix_le_u8 ( const void *  r1,
const void *  r2,
size_t  rsize 
)
static

Definition at line 88 of file mpsort.c.

88  {
89  return _compar_radix_u8(r1, r2, rsize, -1);
90 }

References _compar_radix_u8().

Here is the call graph for this function:

◆ _compar_radix_u8()

static int _compar_radix_u8 ( const void *  r1,
const void *  r2,
size_t  rsize,
int  dir 
)
static

Definition at line 65 of file mpsort.c.

65  {
66  size_t i;
67  /* from most significant */
68  const uint64_t * u1 = (const uint64_t *) r1;
69  const uint64_t * u2 = (const uint64_t *) r2;
70  if(dir < 0) {
71  u1 = (const uint64_t *) ((const char*) u1 + rsize - 8);
72  u2 = (const uint64_t *) ((const char*) u2 + rsize - 8);
73  }
74  for(i = 0; i < rsize; i += 8) {
75  if(*u1 < *u2) return -1;
76  if(*u1 > *u2) return 1;
77  u1 += dir;
78  u2 += dir;
79  }
80  return 0;
81 }

Referenced by _compar_radix_be_u8(), and _compar_radix_le_u8().

Here is the caller graph for this function:

◆ _compute_and_compar_radix()

static int _compute_and_compar_radix ( const void *  p1,
const void *  p2 
)
static

Definition at line 189 of file mpsort.c.

189  {
190  char r1[_cacr_d.rsize], r2[_cacr_d.rsize];
191  _cacr_d.radix(p1, r1, _cacr_d.arg);
192  _cacr_d.radix(p2, r2, _cacr_d.arg);
193  int c1 = _cacr_d.compar(r1, r2, _cacr_d.rsize);
194  return c1;
195 }
static struct crstruct _cacr_d
Definition: mpsort.c:186

References _cacr_d, crstruct::arg, crstruct::compar, crstruct::radix, and crstruct::rsize.

Referenced by radix_sort().

Here is the caller graph for this function:

◆ _create_segment_group()

static void _create_segment_group ( struct SegmentGroupDescr descr,
size_t *  sizes,
size_t  avgsegsize,
int  Ngroup,
MPI_Comm  comm 
)
static

Definition at line 670 of file mpsort.c.

671 {
672  int ThisTask, NTask;
673 
674  MPI_Comm_size(comm, &NTask);
675  MPI_Comm_rank(comm, &ThisTask);
676 
677  descr->ThisSegment = _assign_colors(avgsegsize, sizes, &descr->Nsegments, comm);
678 
679  if(descr->ThisSegment >= 0) {
680  /* assign segments to groups.
681  * if Nsegments < Ngroup, some groups will have no segments, and thus no ranks belong to them. */
682  descr->GroupID = ((size_t) descr->ThisSegment) * Ngroup / descr->Nsegments;
683  } else {
684  descr->GroupID = Ngroup + 1;
685  descr->ThisSegment = NTask + 1;
686  }
687 
688  descr->Ngroup = Ngroup;
689 
690  MPI_Comm_split(comm, descr->GroupID, ThisTask, &descr->Group);
691 
692  MPI_Allreduce(&descr->ThisSegment, &descr->segment_start, 1, MPI_INT, MPI_MIN, descr->Group);
693  MPI_Allreduce(&descr->ThisSegment, &descr->segment_end, 1, MPI_INT, MPI_MAX, descr->Group);
694 
695  descr->segment_end ++;
696 
697  int rank;
698 
699  MPI_Comm_rank(descr->Group, &rank);
700 
701  /* rank with most data in a group is the leader of the group. */
702  descr->group_leader_rank = MPIU_GetLoc(&sizes[ThisTask], MPI_LONG, MPI_MAX, descr->Group);
703 
704  descr->is_group_leader = rank == descr->group_leader_rank;
705 
706  MPI_Comm_split(comm, (rank == descr->group_leader_rank)? 0 : 1, ThisTask, &descr->Leaders);
707 
708  MPI_Comm_split(descr->Group, descr->ThisSegment, ThisTask, &descr->Segment);
709 
710  /* rank with least data in a segment is the leader of the segment. */
711  descr->segment_leader_rank = MPIU_GetLoc(&sizes[ThisTask], MPI_LONG, MPI_MIN, descr->Segment);
712 }
static int MPIU_GetLoc(const void *base, MPI_Datatype type, MPI_Op op, MPI_Comm comm)
Definition: mpsort.c:641
static int _assign_colors(size_t glocalsize, size_t *sizes, int *ncolor, MPI_Comm comm)
Definition: mpsort.c:555
int segment_leader_rank
Definition: mpsort.c:630
MPI_Comm Group
Definition: mpsort.c:631
MPI_Comm Leaders
Definition: mpsort.c:632
int group_leader_rank
Definition: mpsort.c:629
MPI_Comm Segment
Definition: mpsort.c:633
int is_group_leader
Definition: mpsort.c:628

References _assign_colors(), SegmentGroupDescr::Group, SegmentGroupDescr::group_leader_rank, SegmentGroupDescr::GroupID, SegmentGroupDescr::is_group_leader, SegmentGroupDescr::Leaders, MPIU_GetLoc(), SegmentGroupDescr::Ngroup, SegmentGroupDescr::Nsegments, NTask, SegmentGroupDescr::Segment, SegmentGroupDescr::segment_end, SegmentGroupDescr::segment_leader_rank, SegmentGroupDescr::segment_start, SegmentGroupDescr::ThisSegment, and ThisTask.

Referenced by mpsort_mpi_newarray_impl().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _destroy_mpsort_mpi()

static void _destroy_mpsort_mpi ( struct crmpistruct o)
static

Definition at line 523 of file mpsort.c.

523  {
524  MPI_Type_free(&o->MPI_TYPE_RADIX);
525  MPI_Type_free(&o->MPI_TYPE_DATA);
526 }
MPI_Datatype MPI_TYPE_DATA
Definition: mpsort.c:480
MPI_Datatype MPI_TYPE_RADIX
Definition: mpsort.c:479

References crmpistruct::MPI_TYPE_DATA, and crmpistruct::MPI_TYPE_RADIX.

Referenced by mpsort_mpi_newarray_impl().

Here is the caller graph for this function:

◆ _destroy_segment_group()

static void _destroy_segment_group ( struct SegmentGroupDescr descr)
static

Definition at line 715 of file mpsort.c.

716 {
717  MPI_Comm_free(&descr->Segment);
718  MPI_Comm_free(&descr->Group);
719  MPI_Comm_free(&descr->Leaders);
720 }

References SegmentGroupDescr::Group, SegmentGroupDescr::Leaders, and SegmentGroupDescr::Segment.

Referenced by mpsort_mpi_newarray_impl().

Here is the caller graph for this function:

◆ _find_Pmax_Pmin_C()

static void _find_Pmax_Pmin_C ( void *  mybase,
size_t  mynmemb,
size_t  myoutnmemb,
char *  Pmax,
char *  Pmin,
ptrdiff_t *  C,
struct crstruct d,
struct crmpistruct o 
)
static

Definition at line 1304 of file mpsort.c.

1309  {
1310  memset(Pmax, 0, d->rsize);
1311  memset(Pmin, -1, d->rsize);
1312 
1313  char myPmax[d->rsize];
1314  char myPmin[d->rsize];
1315 
1316  size_t * eachnmemb = ta_malloc("eachnmemb", size_t, o->NTask);
1317  size_t * eachoutnmemb = ta_malloc("eachoutnmemb", size_t, o->NTask);
1318  char * eachPmax = (char *) mymalloc("eachPmax", d->rsize * o->NTask * sizeof(char));
1319  char * eachPmin = (char *) mymalloc("eachPmin", d->rsize * o->NTask * sizeof(char));
1320  int i;
1321 
1322  if(mynmemb > 0) {
1323  d->radix((char*) mybase + (mynmemb - 1) * d->size, myPmax, d->arg);
1324  d->radix(mybase, myPmin, d->arg);
1325  } else {
1326  memset(myPmin, 0, d->rsize);
1327  memset(myPmax, 0, d->rsize);
1328  }
1329 
1330  MPI_Allgather(&mynmemb, 1, MPI_TYPE_PTRDIFF,
1331  eachnmemb, 1, MPI_TYPE_PTRDIFF, o->comm);
1332  MPI_Allgather(&myoutnmemb, 1, MPI_TYPE_PTRDIFF,
1333  eachoutnmemb, 1, MPI_TYPE_PTRDIFF, o->comm);
1334  MPI_Allgather(myPmax, 1, o->MPI_TYPE_RADIX,
1335  eachPmax, 1, o->MPI_TYPE_RADIX, o->comm);
1336  MPI_Allgather(myPmin, 1, o->MPI_TYPE_RADIX,
1337  eachPmin, 1, o->MPI_TYPE_RADIX, o->comm);
1338 
1339 
1340  C[0] = 0;
1341  for(i = 0; i < o->NTask; i ++) {
1342  C[i + 1] = C[i] + eachoutnmemb[i];
1343  if(eachnmemb[i] == 0) continue;
1344 
1345  if(d->compar(eachPmax + i * d->rsize, Pmax, d->rsize) > 0) {
1346  memcpy(Pmax, eachPmax + i * d->rsize, d->rsize);
1347  }
1348  if(d->compar(eachPmin + i * d->rsize, Pmin, d->rsize) < 0) {
1349  memcpy(Pmin, eachPmin + i * d->rsize, d->rsize);
1350  }
1351  }
1352 
1353  myfree(eachPmin);
1354  myfree(eachPmax);
1355  myfree(eachoutnmemb);
1356  myfree(eachnmemb);
1357 }
#define mymalloc(name, size)
Definition: mymalloc.h:15
#define ta_malloc(name, type, nele)
Definition: mymalloc.h:25
#define myfree(x)
Definition: mymalloc.h:19
MPI_Comm comm
Definition: mpsort.c:481
int NTask
Definition: mpsort.c:488

References crstruct::arg, crmpistruct::comm, crstruct::compar, piter::d, MPI_TYPE_PTRDIFF, crmpistruct::MPI_TYPE_RADIX, myfree, mymalloc, crmpistruct::NTask, crstruct::radix, crstruct::rsize, crstruct::size, and ta_malloc.

Referenced by mpsort_mpi_histogram_sort().

Here is the caller graph for this function:

◆ _histogram()

static void _histogram ( char *  P,
int  Plength,
void *  mybase,
size_t  mynmemb,
ptrdiff_t *  myCLT,
ptrdiff_t *  myCLE,
struct crstruct d 
)
static

Definition at line 305 of file mpsort.c.

307  {
308  int it;
309 
310  if(myCLT) {
311  myCLT[0] = 0;
312  for(it = 0; it < Plength; it ++) {
313  /* No need to start from the beginging of mybase, since myubase and P are both sorted */
314  ptrdiff_t offset = myCLT[it];
315  myCLT[it + 1] = _bsearch_last_lt(P + it * d->rsize,
316  ((char*) mybase) + offset * d->size,
317  mynmemb - offset, d)
318  + 1 + offset;
319  }
320  myCLT[it + 1] = mynmemb;
321  }
322  if(myCLE) {
323  myCLE[0] = 0;
324  for(it = 0; it < Plength; it ++) {
325  /* No need to start from the beginging of mybase, since myubase and P are both sorted */
326  ptrdiff_t offset = myCLE[it];
327  myCLE[it + 1] = _bsearch_last_le(P + it * d->rsize,
328  ((char*) mybase) + offset * d->size,
329  mynmemb - offset, d)
330  + 1 + offset;
331  }
332  myCLE[it + 1] = mynmemb;
333  }
334 }
static ptrdiff_t _bsearch_last_le(void *P, void *base, size_t nmemb, struct crstruct *d)
Definition: mpsort.c:257
static ptrdiff_t _bsearch_last_lt(void *P, void *base, size_t nmemb, struct crstruct *d)
Definition: mpsort.c:216

References _bsearch_last_le(), _bsearch_last_lt(), P, crstruct::rsize, and crstruct::size.

Referenced by mpsort_mpi_histogram_sort().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _mpsort_mpi_parse_env()

static void _mpsort_mpi_parse_env ( )
static

Definition at line 1424 of file mpsort.c.

1425 {
1426  static int _mpsort_env_parsed = 0;
1427  if(_mpsort_env_parsed) return;
1428 
1429  _mpsort_env_parsed = 1;
1430  if(getenv("MPSORT_DISABLE_GATHER_SORT"))
1432  if(getenv("MPSORT_REQUIRE_GATHER_SORT "))
1434 }
void mpsort_mpi_set_options(int options)
Definition: mpsort.c:1437
#define MPSORT_REQUIRE_GATHER_SORT
Definition: mpsort.h:6
#define MPSORT_DISABLE_GATHER_SORT
Definition: mpsort.h:5

References MPSORT_DISABLE_GATHER_SORT, mpsort_mpi_set_options(), and MPSORT_REQUIRE_GATHER_SORT.

Referenced by mpsort_mpi_has_options(), mpsort_mpi_set_options(), and mpsort_mpi_unset_options().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _setup_mpsort_mpi()

static void _setup_mpsort_mpi ( struct crmpistruct o,
struct crstruct d,
void *  myoutbase,
size_t  myoutnmemb,
MPI_Comm  comm 
)
static

Definition at line 493 of file mpsort.c.

497 {
498 
499  o->comm = comm;
500 
501  MPI_Comm_size(comm, &o->NTask);
502  MPI_Comm_rank(comm, &o->ThisTask);
503 
504  o->mybase = d->base;
505  o->mynmemb = d->nmemb;
506  o->myoutbase = myoutbase;
507  o->myoutnmemb = myoutnmemb;
508 
509  MPI_Allreduce(&o->mynmemb, &o->nmemb, 1, MPI_TYPE_PTRDIFF, MPI_SUM, comm);
510  MPI_Allreduce(&o->myoutnmemb, &o->outnmemb, 1, MPI_TYPE_PTRDIFF, MPI_SUM, comm);
511 
512  if(o->outnmemb != o->nmemb) {
513  endrun(4, "total number of items in the item does not match the input %ld != %ld\n", o->outnmemb, o->nmemb);
514  }
515 
516  MPI_Type_contiguous(d->rsize, MPI_BYTE, &o->MPI_TYPE_RADIX);
517  MPI_Type_commit(&o->MPI_TYPE_RADIX);
518 
519  MPI_Type_contiguous(d->size, MPI_BYTE, &o->MPI_TYPE_DATA);
520  MPI_Type_commit(&o->MPI_TYPE_DATA);
521 
522 }
void endrun(int where, const char *fmt,...)
Definition: endrun.c:147
size_t outnmemb
Definition: mpsort.c:487
size_t myoutnmemb
Definition: mpsort.c:486
void * mybase
Definition: mpsort.c:482
void * myoutbase
Definition: mpsort.c:483
size_t mynmemb
Definition: mpsort.c:484
size_t nmemb
Definition: mpsort.c:485
int ThisTask
Definition: mpsort.c:489
void * base
Definition: mpsort.c:21
size_t nmemb
Definition: mpsort.c:22

References crstruct::base, crmpistruct::comm, endrun(), crmpistruct::MPI_TYPE_DATA, MPI_TYPE_PTRDIFF, crmpistruct::MPI_TYPE_RADIX, crmpistruct::mybase, crmpistruct::mynmemb, crmpistruct::myoutbase, crmpistruct::myoutnmemb, crstruct::nmemb, crmpistruct::nmemb, crmpistruct::NTask, crmpistruct::outnmemb, crstruct::rsize, crstruct::size, and crmpistruct::ThisTask.

Referenced by mpsort_mpi_newarray_impl().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _setup_radix_sort()

void _setup_radix_sort ( struct crstruct d,
void *  base,
size_t  nmemb,
size_t  size,
void(*)(const void *ptr, void *radix, void *arg)  radix,
size_t  rsize,
void *  arg 
)

Definition at line 129 of file mpsort.c.

136  {
137 
138  /* Cheers stack overflow*/
139  union {
140  uint32_t i;
141  char c[4];
142  } be_detect = {0x01020304};
143  d->base = base;
144  d->nmemb = nmemb;
145  d->rsize = rsize;
146  d->arg = arg;
147  d->radix = radix;
148  d->size = size;
149  switch(rsize) {
150  case 2:
151  d->compar = (_compar_fn_t) _compar_radix_uint16_t;
152  d->bisect = (_bisect_fn_t) _bisect_radix_uint16_t;
153  break;
154  case 4:
155  d->compar = (_compar_fn_t) _compar_radix_uint32_t;
156  d->bisect = (_bisect_fn_t) _bisect_radix_uint32_t;
157  break;
158  case 8:
159  d->compar = (_compar_fn_t) _compar_radix_uint64_t;
160  d->bisect = (_bisect_fn_t) _bisect_radix_uint64_t;
161  break;
162  default:
163  if(be_detect.c[0] != 1) {
164  if(rsize % 8 == 0) {
166  } else{
168  }
170  } else {
171  if(rsize % 8 == 0) {
173  } else{
175  }
177  }
178  }
179 }
void(* _bisect_fn_t)(void *r, const void *r1, const void *r2, size_t rsize)
Definition: mpsort.c:18
static int _compar_radix_be_u8(const void *r1, const void *r2, size_t rsize)
Definition: mpsort.c:91
static int _compar_radix_le_u8(const void *r1, const void *r2, size_t rsize)
Definition: mpsort.c:88
static int _compar_radix_le(const void *r1, const void *r2, size_t rsize)
Definition: mpsort.c:82
static void _bisect_radix_le(void *r, const void *r1, const void *r2, size_t rsize)
Definition: mpsort.c:122
static void _bisect_radix_be(void *r, const void *r1, const void *r2, size_t rsize)
Definition: mpsort.c:125
int(* _compar_fn_t)(const void *r1, const void *r2, size_t rsize)
Definition: mpsort.c:17
static int _compar_radix_be(const void *r1, const void *r2, size_t rsize)
Definition: mpsort.c:85
_bisect_fn_t bisect
Definition: mpsort.c:28

Referenced by mpsort_mpi_newarray_impl(), and radix_sort().

Here is the caller graph for this function:

◆ _solve_for_layout_mpi()

static int _solve_for_layout_mpi ( int  NTask,
ptrdiff_t *  C,
ptrdiff_t *  myT_CLT,
ptrdiff_t *  myT_CLE,
ptrdiff_t *  myT_C,
MPI_Comm  comm 
)
static

Definition at line 1360 of file mpsort.c.

1366  {
1367  int i, j;
1368  int ThisTask;
1369  MPI_Comm_rank(comm, &ThisTask);
1370 
1371  /* first assume we just send according to myT_CLT */
1372  for(i = 0; i < NTask; i ++) {
1373  myT_C[i] = myT_CLT[i];
1374  }
1375 
1376  /* Solve for each receiving task i
1377  *
1378  * this solves for GL_C[..., i + 1], which depends on GL_C[..., i]
1379  *
1380  * and we have GL_C[..., 0] == 0 by definition.
1381  *
1382  * this cannot be done in parallel wrt i because of the dependency.
1383  *
1384  * a solution is guaranteed because GL_CLE and GL_CLT
1385  * brackes the total counts C (we've found it with the
1386  * iterative counting.
1387  *
1388  * */
1389 
1390  ptrdiff_t sure = 0;
1391 
1392  /* how many will I surely receive? */
1393  for(j = 0; j < NTask; j ++) {
1394  ptrdiff_t recvcount = myT_C[j];
1395  sure += recvcount;
1396  }
1397  /* let's see if we have enough */
1398  ptrdiff_t deficit = C[ThisTask + 1] - sure;
1399 
1400  for(j = 0; j < NTask; j ++) {
1401  /* deficit solved */
1402  if(deficit == 0) break;
1403  if(deficit < 0) {
1404  endrun(10, "More items than there should be at j=%d: deficit=%ld\n (C: %ld sure %ld)", j, deficit, C[ThisTask+1], sure);
1405  }
1406  /* how much task j can supply ? */
1407  ptrdiff_t supply = myT_CLE[j] - myT_C[j];
1408  if(supply < 0) {
1409  endrun(10, "Less items than there should be at j=%d: supply=%ld (myTCLE %ld myTC %ld)\n", j, supply, myT_CLE[j], myT_C[j]);
1410  }
1411  if(supply <= deficit) {
1412  myT_C[j] += supply;
1413  deficit -= supply;
1414  } else {
1415  myT_C[j] += deficit;
1416  deficit = 0;
1417  }
1418  }
1419 
1420  return 0;
1421 }

References endrun(), NTask, and ThisTask.

Referenced by mpsort_mpi_histogram_sort().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ checksum()

static uint64_t checksum ( void *  base,
size_t  nbytes,
MPI_Comm  comm 
)
static

Definition at line 796 of file mpsort.c.

797 {
798  uint64_t sum = 0;
799  char * ptr = (char *) base;
800  size_t i;
801  for(i = 0; i < nbytes; i ++) {
802  sum += ptr[i];
803  }
804  MPI_Allreduce(MPI_IN_PLACE, &sum, 1, MPI_LONG, MPI_SUM, comm);
805  return sum;
806 }

Referenced by mpsort_mpi_newarray_impl().

Here is the caller graph for this function:

◆ MPIU_Gather()

static void MPIU_Gather ( MPI_Comm  comm,
int  root,
const void *  sendbuffer,
void *  recvbuffer,
int  nsend,
size_t  elsize,
int *  totalnrecv 
)
static

Definition at line 947 of file mpsort.c.

948 {
949  int NTask;
950  int ThisTask;
951 
952  MPI_Comm_size(comm, &NTask);
953  MPI_Comm_rank(comm, &ThisTask);
954 
955  MPI_Datatype dtype;
956  MPI_Type_contiguous(elsize, MPI_BYTE, &dtype);
957  MPI_Type_commit(&dtype);
958 
959  int * recvcount = ta_malloc("recvcount", int, NTask);
960  int * rdispls = ta_malloc("rdispls", int, NTask+1);
961  int i;
962  MPI_Gather(&nsend, 1, MPI_INT, recvcount, 1, MPI_INT, root, comm);
963 
964  rdispls[0] = 0;
965  for(i = 1; i <= NTask; i ++) {
966  rdispls[i] = rdispls[i - 1] + recvcount[i - 1];
967  }
968 
969  if(ThisTask == root) {
970  if(totalnrecv)
971  *totalnrecv = rdispls[NTask];
972  } else {
973  if(totalnrecv)
974  *totalnrecv = 0;
975  }
976 
977  MPI_Gatherv(sendbuffer, nsend, dtype, recvbuffer, recvcount, rdispls, dtype, root, comm);
978 
979  ta_free(rdispls);
980  ta_free(recvcount);
981  MPI_Type_free(&dtype);
982 }
#define ta_free(p)
Definition: mymalloc.h:28

References crmpistruct::comm, NTask, ta_free, ta_malloc, and ThisTask.

Referenced by mpsort_mpi_newarray_impl().

Here is the caller graph for this function:

◆ MPIU_GetLoc()

static int MPIU_GetLoc ( const void *  base,
MPI_Datatype  type,
MPI_Op  op,
MPI_Comm  comm 
)
static

Definition at line 641 of file mpsort.c.

642 {
643  ptrdiff_t lb;
644  ptrdiff_t elsize;
645  MPI_Type_get_extent(type, &lb, &elsize);
646 
647  void * tmp = malloc(elsize);
648  /* find the result of the reduction. */
649  MPI_Allreduce(base, tmp, 1, type, op, comm);
650 
651  int ThisTask;
652  int NTask;
653  MPI_Comm_size(comm, &NTask);
654  MPI_Comm_rank(comm, &ThisTask);
655  int rank = NTask;
656  int ret = -1;
657  if (memcmp(base, tmp, elsize) == 0) {
658  rank = ThisTask;
659  }
660  /* find the rank that is the same as the reduction result */
661  /* avoid MPI_IN_PLACE, since if we are using this code, we have assumed we are using
662  * a crazy MPI impl...
663  * */
664  MPI_Allreduce(&rank, &ret, 1, MPI_INT, MPI_MIN, comm);
665  free(tmp);
666  return ret;
667 }

References NTask, and ThisTask.

Referenced by _create_segment_group().

Here is the caller graph for this function:

◆ MPIU_Scatter()

static void MPIU_Scatter ( MPI_Comm  comm,
int  root,
const void *  sendbuffer,
void *  recvbuffer,
int  nrecv,
size_t  elsize,
int *  totalnsend 
)
static

Definition at line 985 of file mpsort.c.

986 {
987  int NTask;
988  int ThisTask;
989  MPI_Comm_size(comm, &NTask);
990  MPI_Comm_rank(comm, &ThisTask);
991 
992  MPI_Datatype dtype;
993  MPI_Type_contiguous(elsize, MPI_BYTE, &dtype);
994  MPI_Type_commit(&dtype);
995 
996  int * sendcount = ta_malloc("sendcount", int, NTask);
997  int * sdispls = ta_malloc("sdispls", int, NTask+1);
998  int i;
999 
1000  MPI_Gather(&nrecv, 1, MPI_INT, sendcount, 1, MPI_INT, root, comm);
1001 
1002  sdispls[0] = 0;
1003  for(i = 1; i <= NTask; i ++) {
1004  sdispls[i] = sdispls[i - 1] + sendcount[i - 1];
1005  }
1006 
1007  if(ThisTask == root) {
1008  if(totalnsend)
1009  *totalnsend = sdispls[NTask];
1010  } else {
1011  if(totalnsend)
1012  *totalnsend = 0;
1013  }
1014  MPI_Scatterv(sendbuffer, sendcount, sdispls, dtype, recvbuffer, nrecv, dtype, root, comm);
1015 
1016  ta_free(sdispls);
1017  ta_free(sendcount);
1018  MPI_Type_free(&dtype);
1019 }

References crmpistruct::comm, NTask, ta_free, ta_malloc, and ThisTask.

Referenced by mpsort_mpi_newarray_impl().

Here is the caller graph for this function:

◆ mpsort_free_timers()

void mpsort_free_timers ( void  )

Definition at line 745 of file mpsort.c.

746 {
747  if(!(_TIMERS.tmr)) {
748  myfree(_TIMERS.tmr);
749  _TIMERS.tmr = NULL;
750  _TIMERS.ntimer = 0;
751  _TIMERS.curtmr = 0;
752  }
753 }
static struct TIMERS _TIMERS
int ntimer
Definition: mpsort.c:551
struct TIMER * tmr
Definition: mpsort.c:549
int curtmr
Definition: mpsort.c:550

References _TIMERS, TIMERS::curtmr, myfree, TIMERS::ntimer, and TIMERS::tmr.

◆ mpsort_increment_timer()

static void mpsort_increment_timer ( const char *  name,
int  erase 
)
static

Definition at line 723 of file mpsort.c.

724 {
725  if(!(_TIMERS.tmr))
726  return;
727  struct TIMER * tmr = _TIMERS.tmr+_TIMERS.curtmr;
728  if(erase && _TIMERS.curtmr > 0)
729  _TIMERS.curtmr --;
730  tmr->time = MPI_Wtime();
731  strncpy(tmr->name, name, 19);
732  tmr->name[19] = '\0';
733  _TIMERS.curtmr++;
734 }
const char * name
Definition: densitykernel.c:93
Definition: mpsort.c:543
char name[20]
Definition: mpsort.c:545
double time
Definition: mpsort.c:544

References _TIMERS, TIMERS::curtmr, name, TIMER::name, TIMER::time, and TIMERS::tmr.

Referenced by mpsort_mpi_histogram_sort().

Here is the caller graph for this function:

◆ mpsort_mpi_find_ntimers()

int mpsort_mpi_find_ntimers ( struct TIMERS timers)

Definition at line 770 of file mpsort.c.

770  {
771  return _TIMERS.curtmr;
772 }

References _TIMERS, and TIMERS::curtmr.

◆ mpsort_mpi_has_options()

int mpsort_mpi_has_options ( int  options)

Definition at line 1444 of file mpsort.c.

1445 {
1447  return _mpsort_mpi_options & options;
1448 }
static void _mpsort_mpi_parse_env()
Definition: mpsort.c:1424
static int _mpsort_mpi_options
Definition: mpsort.c:460

References _mpsort_mpi_options, and _mpsort_mpi_parse_env().

Referenced by mpsort_mpi_newarray_impl().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ mpsort_mpi_histogram_sort()

static int mpsort_mpi_histogram_sort ( struct crstruct  d,
struct crmpistruct  o 
)
static

Definition at line 1022 of file mpsort.c.

1023 {
1024  ptrdiff_t * myC = (ptrdiff_t *) mymalloc("myhistC", (o.NTask + 1) * sizeof(ptrdiff_t));
1025 
1026  /* Desired counts*/
1027  ptrdiff_t * C = (ptrdiff_t *) mymalloc("histC", (o.NTask + 1) * sizeof(ptrdiff_t));
1028  /* counts of less than P */
1029  ptrdiff_t * myCLT = (ptrdiff_t *) mymalloc("myhistC", (o.NTask + 1) * sizeof(ptrdiff_t));
1030  ptrdiff_t * CLT = (ptrdiff_t *) mymalloc("histCLT", (o.NTask + 1) * sizeof(ptrdiff_t));
1031  /* counts of less than or equal to P */
1032  ptrdiff_t * myCLE = (ptrdiff_t *) mymalloc("myhistCLE", (o.NTask + 1) * sizeof(ptrdiff_t));
1033  ptrdiff_t * CLE = (ptrdiff_t *) mymalloc("CLE", (o.NTask + 1) * sizeof(ptrdiff_t));
1034 
1035  int iter = 0;
1036  int done = 0;
1037  char * buffer;
1038  int i;
1039 
1040  mpsort_increment_timer("START", 0);
1041 
1042  /* and sort the local array */
1043  radix_sort(d.base, d.nmemb, d.size, d.radix, d.rsize, d.arg);
1044 
1045  MPI_Barrier(o.comm);
1046 
1047  mpsort_increment_timer("FirstSort", 0);
1048 
1049  char * P = ta_malloc("PP", char, d.rsize * (o.NTask - 1));
1050  memset(P, 0, d.rsize * (o.NTask -1));
1051 
1052  char Pmax[d.rsize];
1053  char Pmin[d.rsize];
1054 
1055  _find_Pmax_Pmin_C(o.mybase, o.mynmemb, o.myoutnmemb, Pmax, Pmin, C, &d, &o);
1056 
1057  mpsort_increment_timer("PmaxPmin", 0);
1058 
1059  struct piter pi;
1060 
1061  piter_init(&pi, Pmin, Pmax, o.NTask - 1, &d);
1062 
1063  while(!done) {
1064  iter ++;
1065  piter_bisect(&pi, P);
1066 
1067  _histogram(P, o.NTask - 1, o.mybase, o.mynmemb, myCLT, myCLE, &d);
1068 
1069  MPI_Allreduce(myCLT, CLT, o.NTask + 1,
1070  MPI_TYPE_PTRDIFF, MPI_SUM, o.comm);
1071  MPI_Allreduce(myCLE, CLE, o.NTask + 1,
1072  MPI_TYPE_PTRDIFF, MPI_SUM, o.comm);
1073 
1074  char bisectnum[20];
1075  snprintf(bisectnum, 20, "bisect%04d", iter);
1076  mpsort_increment_timer(bisectnum, iter > 10);
1077 
1078  piter_accept(&pi, P, C, CLT, CLE);
1079 #if 0
1080  {
1081  int k;
1082  for(k = 0; k < o.NTask; k ++) {
1083  MPI_Barrier(o.comm);
1084  int i;
1085  if(o.ThisTask != k) continue;
1086 
1087  printf("P (%d): PMin %d PMax %d P ",
1088  o.ThisTask,
1089  *(int*) Pmin,
1090  *(int*) Pmax
1091  );
1092  for(i = 0; i < o.NTask - 1; i ++) {
1093  printf(" %d ", ((int*) P) [i]);
1094  }
1095  printf("\n");
1096 
1097  printf("C (%d): ", o.ThisTask);
1098  for(i = 0; i < o.NTask + 1; i ++) {
1099  printf("%ld ", C[i]);
1100  }
1101  printf("\n");
1102  printf("CLT (%d): ", o.ThisTask);
1103  for(i = 0; i < o.NTask + 1; i ++) {
1104  printf("%ld ", CLT[i]);
1105  }
1106  printf("\n");
1107  printf("CLE (%d): ", o.ThisTask);
1108  for(i = 0; i < o.NTask + 1; i ++) {
1109  printf("%ld ", CLE[i]);
1110  }
1111  printf("\n");
1112 
1113  }
1114  }
1115 #endif
1116  done = piter_all_done(&pi);
1117  }
1118 
1119  piter_destroy(&pi);
1120 
1121  _histogram(P, o.NTask - 1, o.mybase, o.mynmemb, myCLT, myCLE, &d);
1122 
1123  ta_free(P);
1124 
1125  mpsort_increment_timer("findP", 0);
1126 
1127  ptrdiff_t * myT_C = (ptrdiff_t *) mymalloc("myhistT_C", (o.NTask) * sizeof(ptrdiff_t));
1128  ptrdiff_t * myT_CLT = (ptrdiff_t *) mymalloc("myhistCLT", (o.NTask) * sizeof(ptrdiff_t));
1129  ptrdiff_t * myT_CLE = (ptrdiff_t *) mymalloc("myhistCLE", (o.NTask) * sizeof(ptrdiff_t));
1130 
1131  /* transpose the matrix, could have been done with a new datatype */
1132  /*
1133  MPI_Alltoall(myCLT, 1, MPI_TYPE_PTRDIFF,
1134  myT_CLT, 1, MPI_TYPE_PTRDIFF, o.comm);
1135  */
1136  MPI_Alltoall(myCLT + 1, 1, MPI_TYPE_PTRDIFF,
1137  myT_CLT, 1, MPI_TYPE_PTRDIFF, o.comm);
1138 
1139  /*MPI_Alltoall(myCLE, 1, MPI_TYPE_PTRDIFF,
1140  myT_CLE, 1, MPI_TYPE_PTRDIFF, o.comm); */
1141  MPI_Alltoall(myCLE + 1, 1, MPI_TYPE_PTRDIFF,
1142  myT_CLE, 1, MPI_TYPE_PTRDIFF, o.comm);
1143 
1144  mpsort_increment_timer("LayDistr", 0);
1145 
1146  _solve_for_layout_mpi(o.NTask, C, myT_CLT, myT_CLE, myT_C, o.comm);
1147 
1148  myfree(myT_CLE);
1149  myfree(myT_CLT);
1150 
1151  myC[0] = 0;
1152  MPI_Alltoall(myT_C, 1, MPI_TYPE_PTRDIFF,
1153  myC + 1, 1, MPI_TYPE_PTRDIFF, o.comm);
1154 
1155  myfree(myT_C);
1156 
1157 #if 0
1158  for(i = 0;i < o.NTask; i ++) {
1159  int j;
1160  MPI_Barrier(o.comm);
1161  if(o.ThisTask != i) continue;
1162  for(j = 0; j < o.NTask + 1; j ++) {
1163  printf("%d %d %d, ",
1164  myCLT[j],
1165  myC[j],
1166  myCLE[j]);
1167  }
1168  printf("\n");
1169 
1170  }
1171 #endif
1172 
1173 
1174  /* Desired counts*/
1175  myfree(CLE);
1176  myfree(myCLE);
1177  myfree(CLT);
1178  myfree(myCLT);
1179  myfree(C);
1180 
1181  int * SendCount = ta_malloc("SendCount", int, o.NTask);
1182  int * SendDispl = ta_malloc("SendDispl", int, o.NTask);
1183  int * RecvCount = ta_malloc("RecvCount", int, o.NTask);
1184  int * RecvDispl = ta_malloc("RecvDispl", int, o.NTask);
1185 
1186  mpsort_increment_timer("LaySolve", 0);
1187 
1188  for(i = 0; i < o.NTask; i ++) {
1189  SendCount[i] = myC[i + 1] - myC[i];
1190  }
1191 
1192  MPI_Alltoall(SendCount, 1, MPI_INT,
1193  RecvCount, 1, MPI_INT, o.comm);
1194 
1195  SendDispl[0] = 0;
1196  RecvDispl[0] = 0;
1197  size_t totrecv = RecvCount[0];
1198  for(i = 1; i < o.NTask; i ++) {
1199  SendDispl[i] = SendDispl[i - 1] + SendCount[i - 1];
1200  RecvDispl[i] = RecvDispl[i - 1] + RecvCount[i - 1];
1201  if(SendDispl[i] != myC[i]) {
1202  endrun(7, "SendDispl error\n");
1203  }
1204  totrecv += RecvCount[i];
1205  }
1206  if(totrecv != o.myoutnmemb) {
1207  endrun(8, "totrecv = %td, mismatch with %td\n", totrecv, o.myoutnmemb);
1208  }
1209 #if 0
1210  {
1211  int k;
1212  for(k = 0; k < o.NTask; k ++) {
1213  MPI_Barrier(o.comm);
1214 
1215  if(o.ThisTask != k) continue;
1216 
1217  printf("P (%d): ", o.ThisTask);
1218  for(i = 0; i < o.NTask - 1; i ++) {
1219  printf("%d ", ((int*) P) [i]);
1220  }
1221  printf("\n");
1222 
1223  printf("C (%d): ", o.ThisTask);
1224  for(i = 0; i < o.NTask + 1; i ++) {
1225  printf("%d ", C[i]);
1226  }
1227  printf("\n");
1228  printf("CLT (%d): ", o.ThisTask);
1229  for(i = 0; i < o.NTask + 1; i ++) {
1230  printf("%d ", CLT[i]);
1231  }
1232  printf("\n");
1233  printf("CLE (%d): ", o.ThisTask);
1234  for(i = 0; i < o.NTask + 1; i ++) {
1235  printf("%d ", CLE[i]);
1236  }
1237  printf("\n");
1238 
1239  printf("MyC (%d): ", o.ThisTask);
1240  for(i = 0; i < o.NTask + 1; i ++) {
1241  printf("%d ", myC[i]);
1242  }
1243  printf("\n");
1244  printf("MyCLT (%d): ", o.ThisTask);
1245  for(i = 0; i < o.NTask + 1; i ++) {
1246  printf("%d ", myCLT[i]);
1247  }
1248  printf("\n");
1249 
1250  printf("MyCLE (%d): ", o.ThisTask);
1251  for(i = 0; i < o.NTask + 1; i ++) {
1252  printf("%d ", myCLE[i]);
1253  }
1254  printf("\n");
1255 
1256  printf("Send Count(%d): ", o.ThisTask);
1257  for(i = 0; i < o.NTask; i ++) {
1258  printf("%d ", SendCount[i]);
1259  }
1260  printf("\n");
1261  printf("My data(%d): ", o.ThisTask);
1262  for(i = 0; i < mynmemb; i ++) {
1263  printf("%d ", ((int*) mybase)[i]);
1264  }
1265  printf("\n");
1266  }
1267  }
1268 #endif
1269  if(o.myoutbase == o.mybase)
1270  buffer = (char *) mymalloc("mpsortbuffer", d.size * o.myoutnmemb);
1271  else
1272  buffer = (char *) o.myoutbase;
1273 
1275  o.mybase, SendCount, SendDispl, o.MPI_TYPE_DATA,
1276  buffer, RecvCount, RecvDispl, o.MPI_TYPE_DATA,
1277  o.comm);
1278 
1279  if(o.myoutbase == o.mybase) {
1280  memcpy(o.myoutbase, buffer, o.myoutnmemb * d.size);
1281  myfree(buffer);
1282  }
1283 
1284  myfree(RecvDispl);
1285  myfree(RecvCount);
1286  myfree(SendDispl);
1287  myfree(SendCount);
1288 
1289  myfree(myC);
1290  MPI_Barrier(o.comm);
1291  mpsort_increment_timer("Exchange", 0);
1292 
1294 
1295  MPI_Barrier(o.comm);
1296 
1297  mpsort_increment_timer("SecondSort", 0);
1298 
1299  mpsort_increment_timer("End", 0);
1300 
1301  return 0;
1302 }
static void piter_accept(struct piter *pi, char *P, ptrdiff_t *C, ptrdiff_t *CLT, ptrdiff_t *CLE)
Definition: mpsort.c:434
static void radix_sort(void *base, size_t nmemb, size_t size, void(*radix)(const void *ptr, void *radix, void *arg), size_t rsize, void *arg)
Definition: mpsort.c:197
static void _find_Pmax_Pmin_C(void *mybase, size_t mynmemb, size_t myoutnmemb, char *Pmax, char *Pmin, ptrdiff_t *C, struct crstruct *d, struct crmpistruct *o)
Definition: mpsort.c:1304
static int _solve_for_layout_mpi(int NTask, ptrdiff_t *C, ptrdiff_t *myT_CLT, ptrdiff_t *myT_CLE, ptrdiff_t *myT_C, MPI_Comm comm)
Definition: mpsort.c:1360
static void piter_init(struct piter *pi, char *Pmin, char *Pmax, int Plength, struct crstruct *d)
Definition: mpsort.c:344
static void piter_destroy(struct piter *pi)
Definition: mpsort.c:364
static void _histogram(char *P, int Plength, void *mybase, size_t mynmemb, ptrdiff_t *myCLT, ptrdiff_t *myCLE, struct crstruct *d)
Definition: mpsort.c:305
static void mpsort_increment_timer(const char *name, int erase)
Definition: mpsort.c:723
static void piter_bisect(struct piter *pi, char *P)
Definition: mpsort.c:378
static int piter_all_done(struct piter *pi)
Definition: mpsort.c:409
Definition: mpsort.c:336
struct crstruct * d
Definition: mpsort.c:342
int MPI_Alltoallv_smart(void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype sendtype, void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype recvtype, MPI_Comm comm)
Definition: system.c:278

References _find_Pmax_Pmin_C(), _histogram(), _solve_for_layout_mpi(), crstruct::arg, crstruct::base, crmpistruct::comm, piter::d, endrun(), MPI_Alltoallv_smart(), crmpistruct::MPI_TYPE_DATA, MPI_TYPE_PTRDIFF, mpsort_increment_timer(), crmpistruct::mybase, myfree, mymalloc, crmpistruct::mynmemb, crmpistruct::myoutbase, crmpistruct::myoutnmemb, crstruct::nmemb, crmpistruct::NTask, P, piter_accept(), piter_all_done(), piter_bisect(), piter_destroy(), piter_init(), crstruct::radix, radix_sort(), crstruct::rsize, crstruct::size, ta_free, ta_malloc, and crmpistruct::ThisTask.

Referenced by mpsort_mpi_newarray_impl().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ mpsort_mpi_impl()

void mpsort_mpi_impl ( void *  mybase,
size_t  mynmemb,
size_t  size,
void(*)(const void *ptr, void *radix, void *arg)  radix,
size_t  rsize,
void *  arg,
MPI_Comm  comm,
const int  line,
const char *  file 
)

Definition at line 775 of file mpsort.c.

781 {
782  mpsort_mpi_newarray_impl(mybase, mynmemb,
783  mybase, mynmemb,
784  size, radix, rsize, arg, comm, line, file);
785 }
void mpsort_mpi_newarray_impl(void *mybase, size_t mynmemb, void *myoutbase, size_t myoutnmemb, size_t elsize, void(*radix)(const void *ptr, void *radix, void *arg), size_t rsize, void *arg, MPI_Comm comm, const int line, const char *file)
Definition: mpsort.c:809

References mpsort_mpi_newarray_impl().

Here is the call graph for this function:

◆ mpsort_mpi_newarray_impl()

void mpsort_mpi_newarray_impl ( void *  mybase,
size_t  mynmemb,
void *  myoutbase,
size_t  myoutnmemb,
size_t  elsize,
void(*)(const void *ptr, void *radix, void *arg)  radix,
size_t  rsize,
void *  arg,
MPI_Comm  comm,
const int  line,
const char *  file 
)

Definition at line 809 of file mpsort.c.

818 {
819  if(MPI_TYPE_PTRDIFF == 0) {
820  if(MPI_SUCCESS != MPI_Type_match_size(MPI_TYPECLASS_INTEGER, sizeof(ptrdiff_t), &MPI_TYPE_PTRDIFF))
821  endrun(3, "Ptrdiff size %ld not recognised\n", sizeof(ptrdiff_t));
822  }
823 
824  struct SegmentGroupDescr seggrp[1];
825 
826  uint64_t sum1 = checksum(mybase, elsize * mynmemb, comm);
827 
828  int NTask;
829  int ThisTask;
830  MPI_Comm_size(comm, &NTask);
831  MPI_Comm_rank(comm, &ThisTask);
832 
833  if(elsize > 8 && elsize % 8 != 0) {
834  if(ThisTask == 0) {
835  endrun(12, "MPSort: element size is large (%d) but not aligned to 8 bytes. "
836  "This is known to frequently trigger MPI bugs. "
837  "Caller site: %s:%d\n",
838  elsize, file, line);
839  }
840  }
841  if(rsize > 8 && rsize % 8 != 0) {
842  if(ThisTask == 0) {
843  endrun(12, "MPSort: radix size is large (%d) but not aligned to 8 bytes. "
844  "This is known to frequently trigger MPI bugs. "
845  "Caller site: %s:%d\n",
846  rsize, file, line);
847  }
848  }
849 
850  size_t * sizes = ta_malloc("sizes", size_t, NTask);
851  size_t myoffset;
852  size_t totalsize = _collect_sizes(mynmemb, sizes, &myoffset, comm);
853 
854  size_t avgsegsize = NTask; /* combine very small ranks to segments */
855  if (avgsegsize * elsize > 4 * 1024 * 1024) {
856  /* do not use more than 4MB in a segment */
857  avgsegsize = 4 * 1024 * 1024 / elsize;
858  }
860  message(0, "MPSort: gathering all data to a single rank for sorting due to MPSORT_REQUIRE_GATHER_SORT. "
861  "Total number of items is %ld. Caller site: %s:%d\n",
862  totalsize, file, line);
863  avgsegsize = totalsize;
864  }
865 
867  avgsegsize = 0;
868  message(0, "MPSort: disable gathering data into larger chunks due to MPSORT_DISABLE_GATHER_SORT. "
869  "Caller site: %s:%d\n",
870  file, line);
871  }
872 
873  /* use as many groups as possible (some will be empty) but at most 1 segment per group */
874  _create_segment_group(seggrp, sizes, avgsegsize, NTask, comm);
875 
876  myfree(sizes);
877  /* group comm == seg comm */
878 
879  void * mysegmentbase = NULL;
880  void * myoutsegmentbase = NULL;
881  size_t mysegmentnmemb;
882  size_t myoutsegmentnmemb;
883 
884  int groupsize;
885  int grouprank;
886  MPI_Comm_size(seggrp->Group, &groupsize);
887  MPI_Comm_rank(seggrp->Group, &grouprank);
888 
889  MPI_Allreduce(&mynmemb, &mysegmentnmemb, 1, MPI_TYPE_PTRDIFF, MPI_SUM, seggrp->Group);
890  MPI_Allreduce(&myoutnmemb, &myoutsegmentnmemb, 1, MPI_TYPE_PTRDIFF, MPI_SUM, seggrp->Group);
891 
892  if (groupsize > 1) {
893  if(grouprank == seggrp->group_leader_rank) {
894  mysegmentbase = mymalloc("segmentbase", mysegmentnmemb * elsize);
895  myoutsegmentbase = mymalloc("outsegment", myoutsegmentnmemb * elsize);
896  }
897  MPIU_Gather(seggrp->Group, seggrp->group_leader_rank, mybase, mysegmentbase, mynmemb, elsize, NULL);
898  } else {
899  mysegmentbase = mybase;
900  myoutsegmentbase = myoutbase;
901  }
902 
903  /* only do sorting on the group leaders for each segment */
904  if(seggrp->is_group_leader) {
905 
906  struct crstruct d;
907  struct crmpistruct o;
908 
909  _setup_radix_sort(&d, mysegmentbase, mysegmentnmemb, elsize, radix, rsize, arg);
910 
911  _setup_mpsort_mpi(&o, &d, myoutsegmentbase, myoutsegmentnmemb, seggrp->Leaders);
912 
914 
916  }
917 
918  if(groupsize > 1) {
919  MPIU_Scatter(seggrp->Group, seggrp->group_leader_rank, myoutsegmentbase, myoutbase, myoutnmemb, elsize, NULL);
920  }
921 
922 /* {
923  int ntmr;
924  if(seggrp->is_group_leader)
925  ntmr = (mpsort_mpi_find_ntimers(tmr) + 1);
926 
927  MPI_Bcast(&ntmr, 1, MPI_INT, seggrp->group_leader_rank, seggrp->Group);
928  MPI_Bcast(tmr, sizeof(tmr[0]) * ntmr, MPI_BYTE, seggrp->group_leader_rank, seggrp->Group);
929  }*/
930 
931  if(grouprank == seggrp->group_leader_rank) {
932  if(myoutsegmentbase != myoutbase)
933  myfree(myoutsegmentbase);
934  if(mysegmentbase != mybase)
935  myfree(mysegmentbase);
936  }
937 
938  _destroy_segment_group(seggrp);
939 
940  uint64_t sum2 = checksum(myoutbase, elsize * myoutnmemb, comm);
941  if (sum1 != sum2) {
942  endrun(5, "Data changed after sorting; checksum mismatch.\n");
943  }
944 }
void message(int where, const char *fmt,...)
Definition: endrun.c:175
int mpsort_mpi_has_options(int options)
Definition: mpsort.c:1444
static uint64_t checksum(void *base, size_t nbytes, MPI_Comm comm)
Definition: mpsort.c:796
static size_t _collect_sizes(size_t localsize, size_t *sizes, size_t *myoffset, MPI_Comm comm)
Definition: mpsort.c:593
static void MPIU_Gather(MPI_Comm comm, int root, const void *sendbuffer, void *recvbuffer, int nsend, size_t elsize, int *totalnrecv)
Definition: mpsort.c:947
static void _create_segment_group(struct SegmentGroupDescr *descr, size_t *sizes, size_t avgsegsize, int Ngroup, MPI_Comm comm)
Definition: mpsort.c:670
static void _destroy_segment_group(struct SegmentGroupDescr *descr)
Definition: mpsort.c:715
static void _setup_mpsort_mpi(struct crmpistruct *o, struct crstruct *d, void *myoutbase, size_t myoutnmemb, MPI_Comm comm)
Definition: mpsort.c:493
static int mpsort_mpi_histogram_sort(struct crstruct d, struct crmpistruct o)
Definition: mpsort.c:1022
static void _destroy_mpsort_mpi(struct crmpistruct *o)
Definition: mpsort.c:523
static void MPIU_Scatter(MPI_Comm comm, int root, const void *sendbuffer, void *recvbuffer, int nrecv, size_t elsize, int *totalnsend)
Definition: mpsort.c:985
void _setup_radix_sort(struct crstruct *d, void *base, size_t nmemb, size_t size, void(*radix)(const void *ptr, void *radix, void *arg), size_t rsize, void *arg)
Definition: mpsort.c:129
size_t totalsize
Definition: mpsort.c:624

References _collect_sizes(), _create_segment_group(), _destroy_mpsort_mpi(), _destroy_segment_group(), _setup_mpsort_mpi(), _setup_radix_sort(), checksum(), crmpistruct::comm, endrun(), SegmentGroupDescr::Group, SegmentGroupDescr::group_leader_rank, SegmentGroupDescr::is_group_leader, SegmentGroupDescr::Leaders, message(), MPI_TYPE_PTRDIFF, MPIU_Gather(), MPIU_Scatter(), MPSORT_DISABLE_GATHER_SORT, mpsort_mpi_has_options(), mpsort_mpi_histogram_sort(), MPSORT_REQUIRE_GATHER_SORT, crmpistruct::mybase, myfree, mymalloc, crmpistruct::myoutbase, crmpistruct::myoutnmemb, NTask, ta_malloc, ThisTask, and SegmentGroupDescr::totalsize.

Referenced by mpsort_mpi_impl().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ mpsort_mpi_report_last_run()

void mpsort_mpi_report_last_run ( )

Definition at line 756 of file mpsort.c.

756  {
757  if(!(_TIMERS.tmr))
758  return;
759  double last = _TIMERS.tmr[0].time;
760  int i;
761  for(i = 1; i < _TIMERS.curtmr; i++) {
762  struct TIMER * tmr = &_TIMERS.tmr[i];
763  if(0 == strncmp(tmr->name, "END", 20))
764  break;
765  message(0, "%s: %g\n", tmr->name, tmr->time - last);
766  last = tmr->time;
767  }
768 }

References _TIMERS, TIMERS::curtmr, message(), TIMER::name, TIMER::time, and TIMERS::tmr.

Here is the call graph for this function:

◆ mpsort_mpi_set_options()

void mpsort_mpi_set_options ( int  options)

Definition at line 1437 of file mpsort.c.

1438 {
1440  _mpsort_mpi_options |= options;
1441 }

References _mpsort_mpi_options, and _mpsort_mpi_parse_env().

Referenced by _mpsort_mpi_parse_env(), and do_mpsort_test().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ mpsort_mpi_unset_options()

void mpsort_mpi_unset_options ( int  options)

Definition at line 1451 of file mpsort.c.

1452 {
1454  _mpsort_mpi_options &= ~options;
1455 
1456 }

References _mpsort_mpi_options, and _mpsort_mpi_parse_env().

Referenced by do_mpsort_test().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ mpsort_setup_timers()

void mpsort_setup_timers ( int  ntimers)

Definition at line 736 of file mpsort.c.

737 {
738  if(!(_TIMERS.tmr)) {
739  _TIMERS.tmr = (struct TIMER *) mymalloc2("timers", ntimers * sizeof(struct TIMER));
740  _TIMERS.ntimer = ntimers;
741  _TIMERS.curtmr = 0;
742  }
743 }
#define mymalloc2(name, size)
Definition: mymalloc.h:16

References _TIMERS, TIMERS::curtmr, mymalloc2, TIMERS::ntimer, and TIMERS::tmr.

Referenced by do_mpsort_test().

Here is the caller graph for this function:

◆ piter_accept()

static void piter_accept ( struct piter pi,
char *  P,
ptrdiff_t *  C,
ptrdiff_t *  CLT,
ptrdiff_t *  CLE 
)
static

Definition at line 434 of file mpsort.c.

435  {
436  struct crstruct * d = pi->d;
437  int i;
438 #if 0
439  for(i = 0; i < pi->Plength + 1; i ++) {
440  printf("counts %d LT %ld C %ld LE %ld\n",
441  i, CLT[i], C[i], CLE[i]);
442  }
443 #endif
444  for(i = 0; i < pi->Plength; i ++) {
445  if( CLT[i + 1] < C[i + 1] && C[i + 1] <= CLE[i + 1]) {
446  pi->stable[i] = 1;
447  continue;
448  } else {
449  if(CLT[i + 1] >= C[i + 1]) {
450  /* P[i] is too big */
451  memcpy(&pi->Pright[i * d->rsize], &P[i * d->rsize], d->rsize);
452  } else {
453  /* P[i] is too small */
454  memcpy(&pi->Pleft[i * d->rsize], &P[i * d->rsize], d->rsize);
455  }
456  }
457  }
458 }
int * stable
Definition: mpsort.c:337
char * Pright
Definition: mpsort.c:341
char * Pleft
Definition: mpsort.c:340
int Plength
Definition: mpsort.c:339

References piter::d, P, piter::Pleft, piter::Plength, piter::Pright, crstruct::rsize, and piter::stable.

Referenced by mpsort_mpi_histogram_sort().

Here is the caller graph for this function:

◆ piter_all_done()

static int piter_all_done ( struct piter pi)
static

Definition at line 409 of file mpsort.c.

409  {
410  int i;
411  int done = 1;
412 #if 0
413 #pragma omp single
414  for(i = 0; i < pi->Plength; i ++) {
415  printf("P %d stable %d narrow %d\n",
416  i, pi->stable[i], pi->narrow[i]);
417  }
418 #endif
419  for(i = 0; i < pi->Plength; i ++) {
420  if(!pi->stable[i]) {
421  done = 0;
422  break;
423  }
424  }
425  return done;
426 }
int * narrow
Definition: mpsort.c:338

References piter::narrow, piter::Plength, and piter::stable.

Referenced by mpsort_mpi_histogram_sort().

Here is the caller graph for this function:

◆ piter_bisect()

static void piter_bisect ( struct piter pi,
char *  P 
)
static

Definition at line 378 of file mpsort.c.

378  {
379  struct crstruct * d = pi->d;
380  int i;
381  for(i = 0; i < pi->Plength; i ++) {
382  if(pi->stable[i]) continue;
383  if(pi->narrow[i]) {
384  /* The last iteration, test Pright directly */
385  memcpy(&P[i * d->rsize],
386  &pi->Pright[i * d->rsize],
387  d->rsize);
388  pi->stable[i] = 1;
389  } else {
390  /* ordinary iteration */
391  d->bisect(&P[i * d->rsize],
392  &pi->Pleft[i * d->rsize],
393  &pi->Pright[i * d->rsize], d->rsize);
394  /* in case the bisect can't move P beyond left,
395  * the range is too small, so we set flag narrow,
396  * and next iteration we will directly test Pright */
397  if(d->compar(&P[i * d->rsize],
398  &pi->Pleft[i * d->rsize], d->rsize) <= 0) {
399  pi->narrow[i] = 1;
400  }
401  }
402 #if 0
403  printf("bisect %d %u %u %u\n", i, *(int*) &P[i * d->rsize],
404  *(int*) &pi->Pleft[i * d->rsize],
405  *(int*) &pi->Pright[i * d->rsize]);
406 #endif
407  }
408 }

References crstruct::bisect, crstruct::compar, piter::d, piter::narrow, P, piter::Pleft, piter::Plength, piter::Pright, crstruct::rsize, and piter::stable.

Referenced by mpsort_mpi_histogram_sort().

Here is the caller graph for this function:

◆ piter_destroy()

static void piter_destroy ( struct piter pi)
static

Definition at line 364 of file mpsort.c.

364  {
365  myfree(pi->Pright);
366  myfree(pi->Pleft);
367  myfree(pi->narrow);
368  myfree(pi->stable);
369 }

References myfree, piter::narrow, piter::Pleft, piter::Pright, and piter::stable.

Referenced by mpsort_mpi_histogram_sort().

Here is the caller graph for this function:

◆ piter_init()

static void piter_init ( struct piter pi,
char *  Pmin,
char *  Pmax,
int  Plength,
struct crstruct d 
)
static

Definition at line 344 of file mpsort.c.

346  {
347  pi->stable = ta_malloc("stable", int, Plength);
348  memset(pi->stable, 0, Plength * sizeof(int));
349  pi->narrow = ta_malloc("narrow", int, Plength);
350  memset(pi->narrow, 0, Plength * sizeof(int));
351  pi->d = d;
352  pi->Pleft = ta_malloc("left", char, Plength * d->rsize);
353  memset(pi->Pleft, 0, Plength * d->rsize * sizeof(char));
354  pi->Pright = ta_malloc("right", char, Plength * d->rsize);
355  memset(pi->Pright, 0, Plength * d->rsize * sizeof(char));
356  pi->Plength = Plength;
357 
358  int i;
359  for(i = 0; i < pi->Plength; i ++) {
360  memcpy(&pi->Pleft[i * d->rsize], Pmin, d->rsize);
361  memcpy(&pi->Pright[i * d->rsize], Pmax, d->rsize);
362  }
363 }

References piter::d, piter::narrow, piter::Pleft, piter::Plength, piter::Pright, crstruct::rsize, piter::stable, and ta_malloc.

Referenced by mpsort_mpi_histogram_sort().

Here is the caller graph for this function:

◆ radix_sort()

static void radix_sort ( void *  base,
size_t  nmemb,
size_t  size,
void(*)(const void *ptr, void *radix, void *arg)  radix,
size_t  rsize,
void *  arg 
)
static

Definition at line 197 of file mpsort.c.

200  {
201 
202  memset(&_cacr_d, 0, sizeof(struct crstruct));
203  _setup_radix_sort(&_cacr_d, base, nmemb, size, radix, rsize, arg);
204 
206 }
static int _compute_and_compar_radix(const void *p1, const void *p2)
Definition: mpsort.c:189
#define qsort_openmp
Definition: test_exchange.c:14

References _cacr_d, _compute_and_compar_radix(), _setup_radix_sort(), crstruct::arg, crstruct::base, crstruct::nmemb, qsort_openmp, crstruct::radix, crstruct::rsize, and crstruct::size.

Referenced by mpsort_mpi_histogram_sort().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ _cacr_d

struct crstruct _cacr_d
static

Definition at line 129 of file mpsort.c.

Referenced by _compute_and_compar_radix(), and radix_sort().

◆ _mpsort_mpi_options

int _mpsort_mpi_options = 0
static

◆ _TIMERS

struct TIMERS _TIMERS
static

◆ MPI_TYPE_PTRDIFF

MPI_Datatype MPI_TYPE_PTRDIFF = 0
static