QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
comm_quda.h
Go to the documentation of this file.
1 #pragma once
2 #include <cstdint>
3 
4 #ifdef __cplusplus
5 extern "C" {
6 #endif
7 
8  typedef struct MsgHandle_s MsgHandle;
9  typedef struct Topology_s Topology;
10 
11  /* defined in quda.h; redefining here to avoid circular references */
12  typedef int (*QudaCommsMap)(const int *coords, void *fdata);
13 
14  /* implemented in comm_common.cpp */
15 
16  char *comm_hostname(void);
17  double comm_drand(void);
18  Topology *comm_create_topology(int ndim, const int *dims, QudaCommsMap rank_from_coords, void *map_data);
19  void comm_destroy_topology(Topology *topo);
20  int comm_ndim(const Topology *topo);
21  const int *comm_dims(const Topology *topo);
22  const int *comm_coords(const Topology *topo);
23  const int *comm_coords_from_rank(const Topology *topo, int rank);
24  int comm_rank_from_coords(const Topology *topo, const int *coords);
25  int comm_rank_displaced(const Topology *topo, const int displacement[]);
28 
29  // routines related to direct peer-2-peer access
30  void comm_set_neighbor_ranks(Topology *topo=NULL);
31  int comm_neighbor_rank(int dir, int dim);
32 
38  int comm_dim(int dim);
39 
45  int comm_coord(int dim);
46 
56  MsgHandle *comm_declare_send_relative_(const char *func, const char *file, int line,
57  void *buffer, int dim, int dir, size_t nbytes);
58 
59 #define comm_declare_send_relative(buffer, dim, dir, nbytes) \
60  comm_declare_send_relative_(__func__, __FILE__, __LINE__, buffer, dim, dir, nbytes)
61 
71  MsgHandle *comm_declare_receive_relative_(const char *func, const char *file, int line,
72  void *buffer, int dim, int dir, size_t nbytes);
73 
74 #define comm_declare_receive_relative(buffer, dim, dir, nbytes) \
75  comm_declare_receive_relative_(__func__, __FILE__, __LINE__, buffer, dim, dir, nbytes)
76 
88  MsgHandle *comm_declare_strided_send_relative_(const char *func, const char *file, int line,
89  void *buffer, int dim, int dir,
90  size_t blksize, int nblocks, size_t stride);
91 
92 #define comm_declare_strided_send_relative(buffer, dim, dir, blksize, nblocks, stride) \
93  comm_declare_strided_send_relative_(__func__, __FILE__, __LINE__, buffer, dim, dir, blksize, nblocks, stride)
94 
106  MsgHandle *comm_declare_strided_receive_relative_(const char *func, const char *file, int line,
107  void *buffer, int dim, int dir,
108  size_t blksize, int nblocks, size_t stride);
109 
110 #define comm_declare_strided_receive_relative(buffer, dim, dir, blksize, nblocks, stride) \
111  comm_declare_strided_receive_relative_(__func__, __FILE__, __LINE__, buffer, dim, dir, blksize, nblocks, stride)
112 
113  void comm_finalize(void);
114  void comm_dim_partitioned_set(int dim);
115  int comm_dim_partitioned(int dim);
116 
121  int comm_partitioned();
122 
127 
133  const char *comm_dim_partitioned_string(const int *comm_dim_override = 0);
134 
139  const char* comm_dim_topology_string();
140 
147  const char *comm_config_string();
148 
152  void comm_init(int ndim, const int *dims, QudaCommsMap rank_from_coords, void *map_data);
153 
157  void comm_init_common(int ndim, const int *dims, QudaCommsMap rank_from_coords, void *map_data);
158 
162  int comm_rank(void);
163 
167  int comm_size(void);
168 
172  int comm_gpuid(void);
173 
178 
185  void comm_gather_hostname(char *hostname_recv_buf);
186 
192  void comm_gather_gpuid(int *gpuid_recv_buf);
193 
198  void comm_peer2peer_init(const char *hostname_recv_buf);
199 
206  bool comm_peer2peer_present();
207 
213 
220  bool comm_peer2peer_enabled(int dir, int dim);
221 
227  void comm_enable_peer2peer(bool enable);
228 
236  bool comm_intranode_enabled(int dir, int dim);
237 
244  void comm_enable_intranode(bool enable);
245 
249  bool comm_gdr_enabled();
250 
254  bool comm_gdr_blacklist();
255 
263  MsgHandle *comm_declare_send_displaced(void *buffer, const int displacement[], size_t nbytes);
264 
272  MsgHandle *comm_declare_receive_displaced(void *buffer, const int displacement[], size_t nbytes);
273 
283  void *buffer, const int displacement[], size_t blksize, int nblocks, size_t stride);
284 
293  MsgHandle *comm_declare_strided_receive_displaced(void *buffer, const int displacement[],
294  size_t blksize, int nblocks, size_t stride);
295 
296  void comm_free(MsgHandle *&mh);
297  void comm_start(MsgHandle *mh);
298  void comm_wait(MsgHandle *mh);
299  int comm_query(MsgHandle *mh);
300  void comm_allreduce(double* data);
301  void comm_allreduce_max(double* data);
302  void comm_allreduce_min(double* data);
303  void comm_allreduce_array(double* data, size_t size);
304  void comm_allreduce_max_array(double* data, size_t size);
305  void comm_allreduce_int(int* data);
306  void comm_allreduce_xor(uint64_t *data);
307  void comm_broadcast(void *data, size_t nbytes);
308  void comm_barrier(void);
309  void comm_abort(int status);
310 
311  void reduceMaxDouble(double &);
312  void reduceDouble(double &);
313  void reduceDoubleArray(double *, const int len);
314  int commDim(int);
315  int commCoords(int);
316  int commDimPartitioned(int dir);
317  void commDimPartitionedSet(int dir);
318 
325  bool commGlobalReduction();
326  void commGlobalReductionSet(bool global_reduce);
327 
328  bool commAsyncReduction();
329  void commAsyncReductionSet(bool global_reduce);
330 
331 #ifdef __cplusplus
332 }
333 #endif
const int * comm_coords(const Topology *topo)
int comm_rank(void)
Definition: comm_mpi.cpp:82
MsgHandle * comm_declare_strided_receive_relative_(const char *func, const char *file, int line, void *buffer, int dim, int dir, size_t blksize, int nblocks, size_t stride)
void comm_finalize(void)
int commDimPartitioned(int dir)
void comm_destroy_topology(Topology *topo)
bool commAsyncReduction()
MsgHandle * comm_declare_receive_displaced(void *buffer, const int displacement[], size_t nbytes)
Definition: comm_mpi.cpp:130
Topology * comm_create_topology(int ndim, const int *dims, QudaCommsMap rank_from_coords, void *map_data)
Definition: comm_common.cpp:94
void comm_peer2peer_init(const char *hostname_recv_buf)
void comm_abort(int status)
Definition: comm_mpi.cpp:328
MsgHandle * comm_declare_strided_receive_displaced(void *buffer, const int displacement[], size_t blksize, int nblocks, size_t stride)
Definition: comm_mpi.cpp:182
int comm_dim(int dim)
void commDimPartitionedSet(int dir)
int commCoords(int)
int comm_partitioned()
Loop over comm_dim_partitioned(dim) for all comms dimensions.
static int rank
Definition: comm_mpi.cpp:44
void comm_allreduce_max_array(double *data, size_t size)
Definition: comm_mpi.cpp:296
void comm_enable_intranode(bool enable)
Enable / disable intra-node (non-peer-to-peer) communication.
void reduceDoubleArray(double *, const int len)
int comm_coord(int dim)
void comm_allreduce_array(double *data, size_t size)
Definition: comm_mpi.cpp:272
void comm_set_neighbor_ranks(Topology *topo=NULL)
Topology * comm_default_topology(void)
const int * comm_coords_from_rank(const Topology *topo, int rank)
int comm_gpuid(void)
int comm_rank_from_coords(const Topology *topo, const int *coords)
static int size
Definition: comm_mpi.cpp:45
const char * comm_dim_partitioned_string(const int *comm_dim_override=0)
Return a string that defines the comm partitioning (used as a tuneKey)
void comm_gather_hostname(char *hostname_recv_buf)
Gather all hostnames.
Definition: comm_mpi.cpp:47
void comm_enable_peer2peer(bool enable)
Enable / disable peer-to-peer communication: used for dslash policies that do not presently support p...
void comm_init(int ndim, const int *dims, QudaCommsMap rank_from_coords, void *map_data)
Initialize the communications, implemented in comm_single.cpp, comm_qmp.cpp, and comm_mpi.cpp.
Definition: comm_mpi.cpp:58
static int ndim
Definition: layout_hyper.c:53
MsgHandle * comm_declare_strided_send_displaced(void *buffer, const int displacement[], size_t blksize, int nblocks, size_t stride)
Definition: comm_mpi.cpp:153
int comm_size(void)
Definition: comm_mpi.cpp:88
char * comm_hostname(void)
Definition: comm_common.cpp:58
void comm_allreduce_min(double *data)
Definition: comm_mpi.cpp:265
MsgHandle * comm_declare_send_displaced(void *buffer, const int displacement[], size_t nbytes)
Definition: comm_mpi.cpp:107
const char * comm_dim_topology_string()
Return a string that defines the comm topology (for use as a tuneKey)
void comm_start(MsgHandle *mh)
Definition: comm_mpi.cpp:216
int commDim(int)
bool comm_intranode_enabled(int dir, int dim)
void comm_dim_partitioned_set(int dim)
int comm_rank_displaced(const Topology *topo, const int displacement[])
double comm_drand(void)
Definition: comm_common.cpp:82
bool commGlobalReduction()
void comm_free(MsgHandle *&mh)
Definition: comm_mpi.cpp:207
int(* coords)[QUDA_MAX_DIM]
Definition: comm_common.cpp:12
int(* QudaCommsMap)(const int *coords, void *fdata)
Definition: comm_quda.h:12
void comm_init_common(int ndim, const int *dims, QudaCommsMap rank_from_coords, void *map_data)
Initialize the communications common to all communications abstractions.
MsgHandle * comm_declare_strided_send_relative_(const char *func, const char *file, int line, void *buffer, int dim, int dir, size_t blksize, int nblocks, size_t stride)
bool comm_deterministic_reduce()
void commAsyncReductionSet(bool global_reduce)
static int dims[4]
Definition: face_gauge.cpp:41
bool comm_peer2peer_enabled(int dir, int dim)
void comm_allreduce_xor(uint64_t *data)
Definition: comm_mpi.cpp:311
const int * comm_dims(const Topology *topo)
void comm_gather_gpuid(int *gpuid_recv_buf)
Gather all GPU ids.
Definition: comm_mpi.cpp:53
void comm_broadcast(void *data, size_t nbytes)
Definition: comm_mpi.cpp:321
void comm_set_default_topology(Topology *topo)
MsgHandle * comm_declare_receive_relative_(const char *func, const char *file, int line, void *buffer, int dim, int dir, size_t nbytes)
void commDimPartitionedReset()
Reset the comm dim partioned array to zero,.
int comm_ndim(const Topology *topo)
int comm_query(MsgHandle *mh)
Definition: comm_mpi.cpp:228
void comm_allreduce_int(int *data)
Definition: comm_mpi.cpp:304
bool comm_gdr_enabled()
Query if GPU Direct RDMA communication is enabled (global setting)
void comm_set_tunekey_string()
Create the topology and partition strings that are used in tuneKeys.
void comm_wait(MsgHandle *mh)
Definition: comm_mpi.cpp:222
void comm_allreduce(double *data)
Definition: comm_mpi.cpp:242
void reduceDouble(double &)
int comm_neighbor_rank(int dir, int dim)
void comm_allreduce_max(double *data)
Definition: comm_mpi.cpp:258
bool comm_peer2peer_present()
Returns true if any peer-to-peer capability is present on this system (regardless of whether it has b...
int comm_peer2peer_enabled_global()
bool comm_gdr_blacklist()
Query if GPU Direct RDMA communication is blacklisted for this GPU.
void reduceMaxDouble(double &)
MsgHandle * comm_declare_send_relative_(const char *func, const char *file, int line, void *buffer, int dim, int dir, size_t nbytes)
int comm_dim_partitioned(int dim)
const char * comm_config_string()
Return a string that defines the P2P/GDR environment variable configuration (for use as a tuneKey to ...
void comm_barrier(void)
Definition: comm_mpi.cpp:326
void commGlobalReductionSet(bool global_reduce)