74 for (
size_t i = 0; i < num_entries; ++i)
87 void *restrict send_buffer,
88 size_t send_size_asize,
size_t send_size_entry,
89 int tag,
MPI_Comm comm,
int rank_lim,
90 MPI_Request *restrict requests,
91 const int send_size[rank_lim][send_size_asize])
97 for (
size_t rank = 0; rank < (size_t)rank_lim; ++rank)
98 if (send_size[rank][send_size_entry] > 0) {
100 send_size[rank][send_size_entry],
101 MPI_PACKED, (
int)rank, tag,
102 comm, requests + reqOfs),
105 offset += (size_t)send_size[rank][send_size_entry];
107 return (
struct Xt_xmdd_txstat){ .bytes = offset, .num_msg = reqOfs };
112 const struct dist_dir *dst_dist_dir,
113 struct isect **src_dst_intersections,
116 struct isect *src_dst_intersections_ = *src_dst_intersections
119 *
sizeof(**src_dst_intersections));
120 size_t isect_fill = 0;
122 *restrict entries_dst = dst_dist_dir->
entries;
123 size_t num_entries_src = (size_t)src_dist_dir->
num_entries,
124 num_entries_dst = (
size_t)dst_dist_dir->
num_entries;
125 for (
size_t i = 0; i < num_entries_src; ++i)
126 for (
size_t j = 0; j < num_entries_dst; ++j)
130 entries_src[i].
list, entries_dst[j].
list, config);
132 src_dst_intersections_[isect_fill]
136 .idxlist = intersection };
141 *src_dst_intersections
143 isect_fill *
sizeof (*src_dst_intersections_));
147#if __GNUC__ == 11 && __GNUC_MINOR__ <= 2
151#pragma GCC diagnostic push
152#pragma GCC diagnostic ignored "-Wvla-parameter"
158 size_t num_intersections,
159 const struct isect *restrict src_dst_intersections,
160 bool isect_idxlist_delete,
161 size_t send_size_asize,
size_t send_size_idx,
162 int (*send_size)[send_size_asize],
163 unsigned char *buffer,
size_t buf_size,
size_t *ofs,
MPI_Comm comm)
165 int prev_send_rank = -1;
166 size_t num_send_indices_requests = 0;
167 size_t origin = 1 ^ target, ofs_ = *ofs;
169 for (
size_t i = 0; i < num_intersections; ++i)
172 int send_rank = src_dst_intersections[i].rank[target];
173 num_send_indices_requests += send_rank != prev_send_rank;
176 XT_MPI_SEND_BUF_CONST
int *prank
177 = CAST_MPI_SEND_BUF(src_dst_intersections[i].rank + origin);
178 if (send_rank != prev_send_rank && prev_send_rank != -1) {
179 send_size[prev_send_rank][send_size_idx] = position;
180 ofs_ += (size_t)position;
183 prev_send_rank = send_rank;
184 xt_mpi_call(MPI_Pack(prank, 1, MPI_INT, buffer+ofs_, (
int)(buf_size-ofs_),
185 &position, comm), comm);
188 (
int)(buf_size-ofs_), &position, comm);
190 if (isect_idxlist_delete)
193 if (prev_send_rank != -1)
194 send_size[prev_send_rank][send_size_idx] = position;
196 *ofs = ofs_ + (size_t)position;
197 return num_send_indices_requests;
200#if __GNUC__ == 11 && __GNUC_MINOR__ <= 2
201#pragma GCC diagnostic pop
210 - (((csx)b)->start > ((csx)a)->start);
224 struct Xt_com_list *restrict entries = (*dist_dir_results)->entries;
225 size_t num_isect_agg = 0;
227 size_t i = 0, num_shards = (size_t)(*dist_dir_results)->num_entries;
228 while (i < num_shards) {
229 int rank = entries[i].rank;
231 size_t num_stripes = 0;
238 }
while (j < num_shards && entries[j].
rank ==
rank);
241 size_t stripe_ofs = 0;
243 size_t num_stripes_of_intersection
247 num_stripes-stripe_ofs);
249 stripe_ofs += num_stripes_of_intersection;
251 qsort(stripes, num_stripes,
sizeof (*stripes),
stripe_cmp);
253 entries[num_isect_agg].rank = rank;
256 (*dist_dir_results)->num_entries = (int)num_isect_agg;
258 + (
size_t)num_isect_agg
266 const struct isect *a = a_, *b = b_;
274 const struct isect *a = a_, *b = b_;
284 return a->
rank - b->rank;
303 xt_mpi_call(MPI_Comm_test_inter(comm, &is_inter), comm);
309 MPI_Comm merge_comm, local_intra_comm;
310 MPI_Group local_group;
311 xt_mpi_call(MPI_Comm_group(comm, &local_group), comm);
312 xt_mpi_call(MPI_Intercomm_merge(comm, 0, &merge_comm), comm);
313 xt_mpi_call(MPI_Comm_create(merge_comm, local_group, &local_intra_comm),
321 comm, local_intra_comm, config);
323 xt_mpi_call(MPI_Comm_free(&local_intra_comm), local_intra_comm);
add versions of standard API functions not returning on error
#define xrealloc(ptr, size)
struct Xt_stripe * stripes
struct Xt_com_list entries[]
struct Xt_config_ xt_default_config
struct Xt_config_ * Xt_config
implementation of configuration object
struct Xt_xmap_ * Xt_xmap
struct Xt_idxlist_ * Xt_idxlist
void xt_idxlist_get_index_stripes_keep_buf(Xt_idxlist idxlist, struct Xt_stripe *stripes, size_t num_stripes_alloc)
int xt_idxlist_get_num_index_stripes(Xt_idxlist idxlist)
void xt_idxlist_pack(Xt_idxlist idxlist, void *buffer, int buffer_size, int *position, MPI_Comm comm)
Xt_idxlist xt_idxlist_get_intersection_custom(Xt_idxlist idxlist_src, Xt_idxlist idxlist_dst, Xt_config config)
void xt_idxlist_delete(Xt_idxlist idxlist)
Provide non-public declarations common to all index lists.
#define xt_idxlist_get_num_indices(idxlist)
struct Xt_stripes_alloc xt_idxstripes_alloc(size_t num_stripes)
Xt_idxlist xt_idxstripes_congeal(struct Xt_stripes_alloc stripes_alloc)
#define xt_mpi_call(call, comm)
void xt_mpi_comm_mark_exclusive(MPI_Comm comm)
Xt_xmap xt_xmap_dist_dir_intracomm_custom_new(Xt_idxlist src_idxlist, Xt_idxlist dst_idxlist, MPI_Comm comm, Xt_config config)
int xt_com_list_rank_cmp(const void *a_, const void *b_)
int xt_xmdd_cmp_isect_src_rank(const void *a_, const void *b_)
Xt_xmap xt_xmap_dist_dir_new(Xt_idxlist src_idxlist, Xt_idxlist dst_idxlist, MPI_Comm comm)
Xt_xmap xt_xmap_dist_dir_custom_new(Xt_idxlist src_idxlist, Xt_idxlist dst_idxlist, MPI_Comm comm, Xt_config config)
int xt_xmdd_cmp_isect_dst_rank(const void *a_, const void *b_)
void xt_xmdd_free_dist_dirs(struct dist_dir_pair dist_dirs)
size_t xt_xmap_dist_dir_match_src_dst(const struct dist_dir *src_dist_dir, const struct dist_dir *dst_dist_dir, struct isect **src_dst_intersections, Xt_config config)
void xt_xmap_dist_dir_same_rank_merge(struct dist_dir **dist_dir_results)
static void xt_xmdd_free_dist_dir(struct dist_dir *dist_dir)
size_t xt_xmap_dist_dir_pack_intersections(enum xt_xmdd_direction target, size_t num_intersections, const struct isect *restrict src_dst_intersections, bool isect_idxlist_delete, size_t send_size_asize, size_t send_size_idx, int(*send_size)[send_size_asize], unsigned char *buffer, size_t buf_size, size_t *ofs, MPI_Comm comm)
static int stripe_cmp(const void *a, const void *b)
struct Xt_xmdd_txstat xt_xmap_dist_dir_send_intersections(void *restrict send_buffer, size_t send_size_asize, size_t send_size_entry, int tag, MPI_Comm comm, int rank_lim, MPI_Request *restrict requests, const int send_size[rank_lim][send_size_asize])
Utility functions for creation of distributed directories.
Xt_xmap xt_xmap_dist_dir_intercomm_custom_new(Xt_idxlist src_idxlist, Xt_idxlist dst_idxlist, MPI_Comm inter_comm, MPI_Comm intra_comm, Xt_config config)