HavoqGT
havoqgt::mpi Namespace Reference

Namespaces

 detail
 

Classes

class  bfs_priority_queue
 
class  bfs_queue
 
class  bfs_visitor
 
class  delegate_partitioned_graph
 
class  dest_pair_partitioner
 
class  edge_source_partitioner
 
class  edge_target_partitioner
 
class  el_partitioned_t
 
class  get_local_id
 
class  get_owner_id
 
class  high_edge_partitioner
 
class  local_dest_id
 
class  local_source_id
 
class  LogStep
 
class  mailbox_routed
 
class  mpi_communicator
 
class  oned_blocked_partitioned_t
 
struct  OverflowSendInfo
 
class  owner_dest_id
 
class  owner_sort
 
class  owner_source_id
 
class  pr_queue
 
class  pr_visitor
 
class  source_partitioner
 
class  sssp_queue
 
class  sssp_visitor
 
class  termination_detection
 
class  triangle_count_visitor
 
class  triangle_priority_queue
 
class  visitor_queue
 

Typedefs

typedef struct havoqgt::mpi::OverflowSendInfo OverflowSendInfo
 

Functions

template<typename TGraph , typename LevelData , typename ParentData >
void breadth_first_search (TGraph *g, LevelData &level_data, ParentData &parent_data, typename TGraph::vertex_locator s)
 
template<typename EdgeType >
void gen_preferential_attachment_edge_list (std::vector< EdgeType > &local_edges, uint64_t in_base_seed, uint64_t in_node_scale, uint64_t in_edge_scale, double in_beta, double in_prob_rewire, MPI_Comm in_comm)
 
template<typename SegmentManager >
bool operator== (const typename delegate_partitioned_graph< SegmentManager >::edge_iterator &x, const typename delegate_partitioned_graph< SegmentManager >::edge_iterator &y)
 
template<typename SegmentManager >
bool operator!= (const typename delegate_partitioned_graph< SegmentManager >::edge_iterator &x, const typename delegate_partitioned_graph< SegmentManager >::edge_iterator &y)
 
template<typename SegmentManager >
bool operator== (const typename delegate_partitioned_graph< SegmentManager >::vertex_iterator &x, const typename delegate_partitioned_graph< SegmentManager >::vertex_iterator &y)
 
template<typename SegmentManager >
bool operator!= (const typename delegate_partitioned_graph< SegmentManager >::vertex_iterator &x, const typename delegate_partitioned_graph< SegmentManager >::vertex_iterator &y)
 
MPI_Datatype mpi_typeof (char)
 
MPI_Datatype mpi_typeof (signed short)
 
MPI_Datatype mpi_typeof (unsigned char)
 
MPI_Datatype mpi_typeof (unsigned long long)
 
MPI_Datatype mpi_typeof (signed long long)
 
MPI_Datatype mpi_typeof (double)
 
MPI_Datatype mpi_typeof (long double)
 
MPI_Datatype mpi_typeof (std::pair< int, int >)
 
MPI_Datatype mpi_typeof (std::pair< float, int >)
 
MPI_Datatype mpi_typeof (std::pair< double, int >)
 
MPI_Datatype mpi_typeof (std::pair< long double, int >)
 
MPI_Datatype mpi_typeof (std::pair< short, int >)
 
template<typename T >
MPI_Op mpi_typeof (std::greater< T >)
 
template<typename T >
MPI_Op mpi_typeof (std::less< T >)
 
template<typename T >
MPI_Op mpi_typeof (std::plus< T >)
 
template<typename T >
MPI_Op mpi_typeof (std::multiplies< T >)
 
template<typename T >
MPI_Op mpi_typeof (std::logical_and< T >)
 
template<typename T >
MPI_Op mpi_typeof (std::logical_or< T >)
 
void mpi_yield_barrier (MPI_Comm mpi_comm)
 
template<typename T , typename Op >
mpi_all_reduce (T in_d, Op in_op, MPI_Comm mpi_comm)
 
template<typename Vec , typename Op >
void mpi_all_reduce_inplace (Vec &vec, Op in_op, MPI_Comm mpi_comm)
 
template<typename T , typename Op >
void mpi_all_reduce (std::vector< T > &in_vec, std::vector< T > &out_vec, Op in_op, MPI_Comm mpi_comm)
 
template<typename T >
void mpi_all_to_all (std::vector< T > &in_vec, std::vector< int > &in_sendcnts, std::vector< T > &out_vec, std::vector< int > &out_recvcnts, MPI_Comm mpi_comm)
 
template<typename T , typename Partitioner >
void mpi_all_to_all (std::vector< T > &inout_vec, std::vector< T > &temp_vec, Partitioner &owner, MPI_Comm mpi_comm)
 
template<typename T , typename Partitioner >
void mpi_all_to_all_better (std::vector< T > &in_vec, std::vector< T > &out_vec, Partitioner &owner, MPI_Comm mpi_comm)
 
template<typename T >
void mpi_all_to_all_in_place (std::vector< T > &in_out_vec, size_t count, MPI_Comm mpi_comm)
 
template<typename T >
void mpi_all_gather (T _t, std::vector< T > &out_p_vec, MPI_Comm mpi_comm)
 TODO: Add tests. More...
 
template<typename T >
void mpi_all_gather (std::vector< T > &in_send, std::vector< T > &out_recv_gather, MPI_Comm mpi_comm)
 TODO: Add tests, especially with non mpi types, POD. More...
 
template<typename T >
void mpi_all_to_all (std::vector< std::vector< T > > &in_p_vec, std::vector< std::vector< T > > &out_p_vec, MPI_Comm mpi_comm)
 
template<typename T >
void mpi_bcast (T &data, int root, MPI_Comm comm)
 
std::ostream & cout_rank0 ()
 
std::ostream & cout_rank0_barrier ()
 
int mpi_comm_rank ()
 
int mpi_comm_size ()
 
template<typename TGraph , typename PRData >
void page_rank (TGraph &g, PRData &pr_data)
 
template<typename TGraph , typename PathData , typename EdgeWeight >
void single_source_shortest_path (TGraph &g, PathData &path_data, EdgeWeight &edge_data, typename TGraph::vertex_locator s)
 
template<typename TGraph >
uint64_t triangle_count (TGraph &g, typename TGraph::vertex_locator s)
 
template<typename Container >
void free_edge_container (Container &edges)
 Frees the container of edges. More...
 
template<>
void free_edge_container< std::vector< std::pair< uint64_t, uint64_t > > > (std::vector< std::pair< uint64_t, uint64_t > > &edges)
 

Typedef Documentation

This class is used to determine where to send a high edge. If the edge's destination is owned by another node, then the edge is sent to that node. Otherwise, it is sent to a node based on the transfer_info.

Transfer_info is a map of delgate ids (not vertex ids) to a dequeue. i.e. delgate_id = m_map_delegate_locator(vertex_id) Each dequeue contains one or more OverflowSendInfo object which contains to_send_id, to_send_count, temp_to_send_count. which determines who will recieve the extra edges.

Function Documentation

template<typename TGraph , typename LevelData , typename ParentData >
void havoqgt::mpi::breadth_first_search ( TGraph *  g,
LevelData &  level_data,
ParentData &  parent_data,
typename TGraph::vertex_locator  s 
)

Definition at line 251 of file breadth_first_search.hpp.

254  {
255 
256  typedef bfs_visitor<TGraph, LevelData, ParentData> visitor_type;
257  visitor_type::set_level_data(&level_data);
258  visitor_type::set_parent_data(&parent_data);
259  typedef visitor_queue< visitor_type, bfs_priority_queue, TGraph > visitor_queue_type;
260 
261  visitor_queue_type vq(g);
262  vq.init_visitor_traversal(s);
263 }

Here is the caller graph for this function:

std::ostream& havoqgt::mpi::cout_rank0 ( )
inline

Definition at line 565 of file mpi.hpp.

565  {
566  int mpi_rank(0);
567  CHK_MPI( MPI_Comm_rank( MPI_COMM_WORLD, &mpi_rank) );
568 
569  if(mpi_rank == 0) {
570  return std::cout;
571  }
573 }
std::ostream & get_null_ostream()
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the call graph for this function:

std::ostream& havoqgt::mpi::cout_rank0_barrier ( )
inline

Definition at line 575 of file mpi.hpp.

575  {
576  CHK_MPI( MPI_Barrier( MPI_COMM_WORLD ));
577  int mpi_rank(0);
578  CHK_MPI( MPI_Comm_rank( MPI_COMM_WORLD, &mpi_rank) );
579 
580  if(mpi_rank == 0) {
581  return std::cout;
582  }
584 }
std::ostream & get_null_ostream()
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the call graph for this function:

template<typename Container >
void havoqgt::mpi::free_edge_container ( Container &  edges)

Frees the container of edges.

Definition at line 60 of file utilities.hpp.

60 {};
template<>
void havoqgt::mpi::free_edge_container< std::vector< std::pair< uint64_t, uint64_t > > > ( std::vector< std::pair< uint64_t, uint64_t > > &  edges)

Definition at line 63 of file utilities.hpp.

63  {
64  std::vector< std::pair<uint64_t, uint64_t> >empty(0);
65  edges.swap(empty);
66 };
template<typename EdgeType >
void havoqgt::mpi::gen_preferential_attachment_edge_list ( std::vector< EdgeType > &  local_edges,
uint64_t  in_base_seed,
uint64_t  in_node_scale,
uint64_t  in_edge_scale,
double  in_beta,
double  in_prob_rewire,
MPI_Comm  in_comm 
)

Definition at line 65 of file gen_preferential_attachment_edge_list.hpp.

72 {
73  namespace ad = havoqgt::detail;
74  int mpi_rank, mpi_size;
75  CHK_MPI( MPI_Comm_rank(in_comm, &mpi_rank) );
76  CHK_MPI( MPI_Comm_size(in_comm, &mpi_size) );
77 
78  //
79  // Calc global number of nodes and edges
80  uint64_t global_num_nodes = uint64_t(1) << in_node_scale;
81  uint64_t global_num_edges = uint64_t(1) << in_edge_scale;
82  uint64_t k = global_num_edges / global_num_nodes;
83  uint64_t edges_per_rank = global_num_edges / mpi_size;
84 
85  ad::preferential_attachment_helper<uint64_t> pa(k, global_num_edges, in_beta,
86  in_base_seed*mpi_rank + mpi_rank);
87 
88  //
89  // Generate inital pa edges. Includes unresolved links
90  double init_time_start = MPI_Wtime();
91  local_edges.resize(edges_per_rank);
92  for(uint64_t i=0; i<edges_per_rank; ++i) {
93  uint64_t edge_index = uint64_t(mpi_rank) + i*uint64_t(mpi_size);
94  local_edges[i] = pa.gen_edge(edge_index);
95  }
96  double init_time_end = MPI_Wtime();
97  if(mpi_rank == 0) {
98  std::cout << "Initial rng for edges took " << init_time_end - init_time_start
99  << " seconds. " << std::endl;
100  }
101 
102  //
103  // We loop until all edges are resolved
104  uint64_t iteration_count(1), count_missing(0), count_total_mising(0);
105  do {
106  //
107  // Count number still missing
108  for(uint64_t i=0; i<edges_per_rank; ++i) {
109  if(pa.is_pointer(local_edges[i].second)) {
110  ++count_total_mising;
111  }
112  }
113 
114  count_missing = 0;
115  double time_start = MPI_Wtime();
116  //
117  // Generate vectors to exchange
118  std::vector<std::vector<uint64_t> > tmp_send_to_rank(mpi_size);
119  std::vector<std::vector<uint64_t> > tmp_mark_send_to_rank(mpi_size);
120  uint64_t count_to_send(0);
121  for(uint64_t i=0; i<local_edges.size(); ++i) {
122  if(pa.is_pointer(local_edges[i].second)) {
123  int owner = pa.value_of_pointer(local_edges[i].second) % mpi_size;
124  tmp_mark_send_to_rank[owner].push_back(i);
125  tmp_send_to_rank[owner].push_back(pa.value_of_pointer(local_edges[i].second));
126  ++count_to_send;
127  }
128  }
129  //
130  // Combine send/recv vectors
131  std::vector<uint64_t> to_send(count_to_send); to_send.reserve(1);
132  std::vector<uint64_t> mark_to_send(count_to_send);
133  std::vector<int> sendcnts(mpi_size,0); sendcnts.reserve(1);
134  to_send.clear(); mark_to_send.clear();
135  for(int i = 0; i < mpi_size; ++i) {
136  for(size_t j = 0; j < tmp_send_to_rank[i].size(); ++j) {
137  to_send.push_back(tmp_send_to_rank[i][j]);
138  mark_to_send.push_back(tmp_mark_send_to_rank[i][j]);
139  }
140  sendcnts[i] = tmp_send_to_rank[i].size();
141  }
142  //
143  // Release memory from temp vectors
144  {
145  std::vector<std::vector<uint64_t> > empty;
146  tmp_send_to_rank.swap(empty);
147  }
148  {
149  std::vector<std::vector<uint64_t> > empty;
150  tmp_mark_send_to_rank.swap(empty);
151  }
152  //
153  // Exchange vectors
154  std::vector<uint64_t> to_recv; // not sized
155  std::vector<int> recvcnts; // not sized
156  havoqgt::mpi::mpi_all_to_all(to_send, sendcnts, to_recv, recvcnts, MPI_COMM_WORLD);
157  //
158  // Look up pointers, pointer jump!
159  for(size_t i=0; i<to_recv.size(); ++i) {
160  uint64_t local_index = to_recv[i]/mpi_size;
161  to_recv[i] = local_edges[local_index].second;
162  }
163  //
164  // Return exchange vectors
165  havoqgt::mpi::mpi_all_to_all(to_recv, recvcnts, to_send, sendcnts, MPI_COMM_WORLD);
166  //
167  // Return new pointers to marked location
168  for(size_t i=0; i<to_send.size(); ++i) {
169  local_edges[ mark_to_send[i] ].second = to_send[i];;
170  }
171  //
172  // Count number still missing
173  for(uint64_t i=0; i<edges_per_rank; ++i) {
174  if(pa.is_pointer(local_edges[i].second)) {
175  ++count_missing;
176  }
177  }
178  count_missing = havoqgt::mpi::mpi_all_reduce(count_missing, std::plus<uint64_t>(), MPI_COMM_WORLD);
179 
180  //
181  // Iteration ouput
182  double time_end = MPI_Wtime();
183  if(mpi_rank == 0) {
184  std::cout << "Iteration " << iteration_count << " took "
185  << time_end - time_start << " seconds. "
186  << count_missing << " still missing. " << std::endl;
187  }
188  ++iteration_count;
189  } while(count_missing > 0);
190  //
191  // Output total nuumber of global exchanges
192  count_total_mising = havoqgt::mpi::mpi_all_reduce(count_total_mising, std::plus<uint64_t>(), MPI_COMM_WORLD);
193  if(mpi_rank == 0) {
194  std::cout << "Total missing (how much data globally exchanged) = " << count_total_mising << std::endl;
195  }
196 
197 
198  //
199  // Randomly rewire
200  if(in_prob_rewire > double(0)) {
201  boost::mt19937 rng(in_base_seed + uint64_t(mpi_rank) * 3ULL);
202  boost::random::uniform_int_distribution<> rand_node(0,global_num_nodes-1);
203  boost::random::uniform_01<> rand_real;
204  for(size_t i=0; i<local_edges.size(); ++i) {
205  if(rand_real(rng) < in_prob_rewire) {
206  EdgeType rand_edge(rand_node(rng), rand_node(rng));
207  local_edges[i] = rand_edge;
208  }
209  }
210  }
211 
212  //
213  // Scramble edges
214  for(size_t i=0; i<local_edges.size(); ++i) {
215  //TODO: This needs a mod because we never exactly make the correct
216  //number of vertices. We need a better solution!
217  local_edges[i].first %= global_num_nodes;
218  local_edges[i].second %= global_num_nodes;
219  local_edges[i].first = ad::hash_nbits(local_edges[i].first, in_node_scale);
220  local_edges[i].second = ad::hash_nbits(local_edges[i].second, in_node_scale);
221  }
222 
223  // //
224  // // Symmetrizing b/c PA is always undirected.
225  // uint64_t old_size = local_edges.size();
226  // local_edges.reserve(old_size * 2);
227  // for(uint64_t i=0; i<old_size; ++i) {
228  // EdgeType edge;
229  // edge.first = local_edges[i].second;
230  // edge.second = local_edges[i].first;
231  // local_edges.push_back( edge );
232  // }
233 
234  // //
235  // // Shuffle edges because we need to mix the hubs in for future partitioning
236  // std::random_shuffle(local_edges.begin(), local_edges.end());
237 }
uint64_t hash_nbits(uint64_t input, int n)
Definition: hash.hpp:115
T mpi_all_reduce(T in_d, Op in_op, MPI_Comm mpi_comm)
Definition: mpi.hpp:176
#define CHK_MPI(a)
Definition: mpi.hpp:68
void mpi_all_to_all(std::vector< T > &in_vec, std::vector< int > &in_sendcnts, std::vector< T > &out_vec, std::vector< int > &out_recvcnts, MPI_Comm mpi_comm)
Definition: mpi.hpp:218

Here is the call graph for this function:

template<typename T >
void havoqgt::mpi::mpi_all_gather ( _t,
std::vector< T > &  out_p_vec,
MPI_Comm  mpi_comm 
)

TODO: Add tests.

Definition at line 432 of file mpi.hpp.

432  {
433  int mpi_size(0);
434  CHK_MPI( MPI_Comm_size( mpi_comm, &mpi_size) );
435 
436  out_p_vec.resize(mpi_size);
437  CHK_MPI( MPI_Allgather(&_t, sizeof(_t), MPI_BYTE,
438  &(out_p_vec[0]), sizeof(_t), MPI_BYTE, mpi_comm) );
439 }
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the caller graph for this function:

template<typename T >
void havoqgt::mpi::mpi_all_gather ( std::vector< T > &  in_send,
std::vector< T > &  out_recv_gather,
MPI_Comm  mpi_comm 
)

TODO: Add tests, especially with non mpi types, POD.

Definition at line 443 of file mpi.hpp.

444  {
445 
446  int mpi_size(0), mpi_rank(0);
447  CHK_MPI( MPI_Comm_size( mpi_comm, &mpi_size) );
448  CHK_MPI( MPI_Comm_rank( mpi_comm, &mpi_rank) );
449 
450  int my_size = in_send.size();
451  std::vector<int> recv_counts(mpi_size,0);
452  std::vector<int> recv_disps(mpi_size,0);
453 
454  //Gather recv counts for all ranks
455  mpi_all_gather(my_size, recv_counts, mpi_comm);
456 
457  //Allocate recv vector
458  int total_recv = std::accumulate(recv_counts.begin(), recv_counts.end(), 0);
459  if(total_recv > 0) {
460  out_recv_gather.resize(total_recv);
461 
462  //cacl recv disps
463  std::partial_sum(recv_counts.begin(), recv_counts.end(), recv_disps.begin());
464  for(size_t i=0; i<recv_disps.size(); ++i) {
465  recv_disps[i] -= recv_counts[i]; //set to 0 offset
466  }
467 
468  void* send_buff = in_send.size() == 0 ? NULL : &(in_send[0]);
469  CHK_MPI( MPI_Allgatherv( send_buff, my_size, mpi_typeof(T()),
470  &(out_recv_gather[0]), &(recv_counts[0]), &(recv_disps[0]),
471  mpi_typeof(T()), mpi_comm));
472  }
473 }
void mpi_all_gather(std::vector< T > &in_send, std::vector< T > &out_recv_gather, MPI_Comm mpi_comm)
TODO: Add tests, especially with non mpi types, POD.
Definition: mpi.hpp:443
MPI_Op mpi_typeof(std::logical_or< T >)
Definition: mpi.hpp:114
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the call graph for this function:

template<typename T , typename Op >
T havoqgt::mpi::mpi_all_reduce ( in_d,
Op  in_op,
MPI_Comm  mpi_comm 
)

Definition at line 176 of file mpi.hpp.

176  {
177  T to_return;
178  CHK_MPI(
179  MPI_Allreduce( &in_d, &to_return, 1, mpi_typeof(T()), mpi_typeof(in_op),
180  mpi_comm)
181  );
182  return to_return;
183 }
MPI_Op mpi_typeof(std::logical_or< T >)
Definition: mpi.hpp:114
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the call graph for this function:

Here is the caller graph for this function:

template<typename T , typename Op >
void havoqgt::mpi::mpi_all_reduce ( std::vector< T > &  in_vec,
std::vector< T > &  out_vec,
Op  in_op,
MPI_Comm  mpi_comm 
)

Definition at line 194 of file mpi.hpp.

195  {
196  out_vec.resize(in_vec.size());
197  CHK_MPI(
198  MPI_Allreduce( &(in_vec[0]), &(out_vec[0]), in_vec.size(),
199  mpi_typeof(in_vec[0]), mpi_typeof(in_op), mpi_comm)
200  );
201 }
MPI_Op mpi_typeof(std::logical_or< T >)
Definition: mpi.hpp:114
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the call graph for this function:

template<typename Vec , typename Op >
void havoqgt::mpi::mpi_all_reduce_inplace ( Vec &  vec,
Op  in_op,
MPI_Comm  mpi_comm 
)

Definition at line 186 of file mpi.hpp.

186  {
187  CHK_MPI(
188  MPI_Allreduce(MPI_IN_PLACE, &(vec[0]), vec.size(),
189  mpi_typeof(typename Vec::value_type()), mpi_typeof(in_op), mpi_comm)
190  );
191 }
MPI_Op mpi_typeof(std::logical_or< T >)
Definition: mpi.hpp:114
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the call graph for this function:

Here is the caller graph for this function:

template<typename T >
void havoqgt::mpi::mpi_all_to_all ( std::vector< T > &  in_vec,
std::vector< int > &  in_sendcnts,
std::vector< T > &  out_vec,
std::vector< int > &  out_recvcnts,
MPI_Comm  mpi_comm 
)

Definition at line 218 of file mpi.hpp.

220  {
221  int mpi_size(0), mpi_rank(0);
222  CHK_MPI( MPI_Comm_size( mpi_comm, &mpi_size) );
223  CHK_MPI( MPI_Comm_rank( mpi_comm, &mpi_rank) );
224 
225  void* send_vec = in_vec.size() > 0 ? (void*) &in_vec[0] : NULL;
226  int* send_cnts = in_sendcnts.size() > 0 ? &in_sendcnts[0] : NULL;
227 
228  //
229  // Exchange send counts
230  out_recvcnts.resize(mpi_size);
231  CHK_MPI( MPI_Alltoall( (void*) send_cnts, sizeof(int), MPI_BYTE,
232  (void*) &(out_recvcnts[0]), sizeof(int), MPI_BYTE,
233  mpi_comm ) );
234 
235  //
236  // Calc send & recv disps
237  std::vector<int> sdispls(mpi_size,0), rdispls(mpi_size,0);
238  std::partial_sum(in_sendcnts.begin(), in_sendcnts.end(), sdispls.begin());
239  for(size_t i=0; i<sdispls.size(); ++i) {
240  sdispls[i] -= in_sendcnts[i]; //set to 0 offset
241  }
242  std::partial_sum(out_recvcnts.begin(), out_recvcnts.end(), rdispls.begin());
243  for(size_t i=0; i<rdispls.size(); ++i) {
244  rdispls[i] -= out_recvcnts[i]; //set to 0 offset
245  }
246 
247  out_vec.resize(std::accumulate(out_recvcnts.begin(), out_recvcnts.end(), 0));
248 
249  int* send_displs = sdispls.size() > 0 ? &sdispls[0] : NULL;
250  CHK_MPI( MPI_Alltoallv(send_vec, send_cnts, send_displs,
251  mpi_typeof(T()), &(out_vec[0]), &(out_recvcnts[0]),
252  &(rdispls[0]), mpi_typeof(T()), mpi_comm) );
253 
254 }
MPI_Op mpi_typeof(std::logical_or< T >)
Definition: mpi.hpp:114
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the call graph for this function:

Here is the caller graph for this function:

template<typename T , typename Partitioner >
void havoqgt::mpi::mpi_all_to_all ( std::vector< T > &  inout_vec,
std::vector< T > &  temp_vec,
Partitioner &  owner,
MPI_Comm  mpi_comm 
)

Definition at line 258 of file mpi.hpp.

259  {
260  int mpi_size(0), mpi_rank(0);
261  CHK_MPI( MPI_Comm_size( mpi_comm, &mpi_size) );
262  CHK_MPI( MPI_Comm_rank( mpi_comm, &mpi_rank) );
263 
264  std::vector<int> send_counts(mpi_size,0);
265  std::vector<int> send_disps(mpi_size,0);
266 
267  std::vector<int> recv_counts(mpi_size,0);
268  std::vector<int> recv_disps(mpi_size,0);
269 
270  //sort send vector by owner
271  //std::sort(inout_vec.begin(), inout_vec.end(), owner_sort<T,Partitioner>(owner));
272 
273  //calc send counts
274  for(size_t i=0; i<inout_vec.size(); ++i) {
275  send_counts[owner(inout_vec[i])] += sizeof(T); //sizeof(t) lets us use PODS
276  //if(i>0) {
277  // assert(owner(inout_vec[i-1]) <= owner(inout_vec[i]));
278  //}
279  }
280 
281  //cacl send disps
282  std::partial_sum(send_counts.begin(), send_counts.end(), send_disps.begin());
283  for(size_t i=0; i<send_disps.size(); ++i) {
284  send_disps[i] -= send_counts[i]; //set to 0 offset
285  }
286 
287  temp_vec.resize(inout_vec.size());
288  for(size_t i=0; i<inout_vec.size(); ++i) {
289  std::vector<int> temp_arrange(mpi_size,0);
290  int dest_rank = owner(inout_vec[i]);
291  size_t dest_offset = send_disps[dest_rank] + temp_arrange[dest_rank];
292  temp_arrange[dest_rank] += sizeof(T);
293  dest_offset /= sizeof(T);
294  temp_vec[dest_offset] = inout_vec[i];
295  }
296 
297 
298  //exchange send counts
299  CHK_MPI( MPI_Alltoall( (void*) &(send_counts[0]), sizeof(int), MPI_BYTE,
300  (void*) &(recv_counts[0]), sizeof(int), MPI_BYTE,
301  mpi_comm ) );
302 
303  //Allocate recv vector
304  int total_recv = std::accumulate(recv_counts.begin(), recv_counts.end(), 0);
305  inout_vec.resize(total_recv/sizeof(T));
306 
307  //cacl recv disps
308  std::partial_sum(recv_counts.begin(), recv_counts.end(), recv_disps.begin());
309  for(size_t i=0; i<recv_disps.size(); ++i) {
310  recv_disps[i] -= recv_counts[i]; //set to 0 offset
311  }
312 
313  //perform actual a3lloallv
314  CHK_MPI( MPI_Alltoallv( (void*) &(temp_vec[0]), &(send_counts[0]), &(send_disps[0]), MPI_BYTE,
315  (void*) &(inout_vec[0]),&(recv_counts[0]), &(recv_disps[0]), MPI_BYTE,
316  mpi_comm) );
317 }
#define CHK_MPI(a)
Definition: mpi.hpp:68
template<typename T >
void havoqgt::mpi::mpi_all_to_all ( std::vector< std::vector< T > > &  in_p_vec,
std::vector< std::vector< T > > &  out_p_vec,
MPI_Comm  mpi_comm 
)

All to All exchange of std::vector<T>

Parameters
[in]in_p_vecinput vector
[out]out_p_vecoutput vector
[in]mpi_commMPI communicator
Todo:
Add ability to optimize for sparse sends. Put MPI_Request objects into a resizing vecor to be skipped for 0 byte send/recvs

Definition at line 487 of file mpi.hpp.

489  {
490  int mpi_size(0), mpi_rank(0);
491  CHK_MPI( MPI_Comm_size( mpi_comm, &mpi_size) );
492  CHK_MPI( MPI_Comm_rank( mpi_comm, &mpi_rank) );
493  assert( mpi_size == (int) in_p_vec.size() );
494 
495  std::vector<size_t> per_rank_send_counts(mpi_size);
496  std::vector<size_t> per_rank_recv_counts(mpi_size);
497 
498  for(int i=0; i<mpi_size; ++i) {
499  per_rank_send_counts[i] = in_p_vec[i].size();
500  }
501 
502  CHK_MPI( MPI_Alltoall( &(per_rank_send_counts[0]), 1, mpi_typeof(size_t()),
503  &(per_rank_recv_counts[0]), 1, mpi_typeof(size_t()),
504  mpi_comm ) );
505  //Allocate recv memory
506  out_p_vec.resize(mpi_size);
507  for(int i=0; i<mpi_size; ++i) {
508  out_p_vec[i].resize(per_rank_recv_counts[i]);
509  }
510 
511 /*
512  * //Aggregisive method, good for small-med p?
513  * std::vector<MPI_Request> vec_req;
514  * for(int i=0; i<mpi_size; ++i) {
515  * int send_to_rank = (mpi_rank + i) % mpi_size;
516  * int recv_from_rank = (mpi_rank - i + mpi_size) % mpi_size;
517  * //Post Irecvs
518  * if(out_p_vec[recv_from_rank].size() > 0) {
519  * MPI_Request req;
520  * CHK_MPI( MPI_Irecv( &(out_p_vec[recv_from_rank][0]),
521  * (out_p_vec[recv_from_rank].size() * sizeof(T)),
522  * MPI_BYTE,
523  * recv_from_rank, 0, mpi_comm, &req ) );
524  * vec_req.push_back(req);
525  * }
526  *
527  * //Post Isends
528  * if(in_p_vec[send_to_rank].size() > 0) {
529  * MPI_Request req;
530  * CHK_MPI( MPI_Isend( &(in_p_vec[send_to_rank][0]),
531  * (in_p_vec[send_to_rank].size() * sizeof(T)),
532  * MPI_BYTE,
533  * send_to_rank, 0, mpi_comm, &req) );
534  * vec_req.push_back(req);
535  * }
536  * }
537  * CHK_MPI( MPI_Waitall(vec_req.size(), &(vec_req[0]), MPI_STATUSES_IGNORE) );
538  */
539 
540  //Basic method -- good for large p?
541  //For each rank, in parallel do:
542  for(int i=0; i<mpi_size; ++i) {
543  MPI_Request request[2];
544  int send_to_rank = (mpi_rank + i) % mpi_size;
545  int recv_from_rank = (mpi_rank - i + mpi_size) % mpi_size;
546  CHK_MPI( MPI_Isend( &(in_p_vec[send_to_rank][0]),
547  (in_p_vec[send_to_rank].size() * sizeof(T)),
548  MPI_BYTE,
549  send_to_rank, 0, mpi_comm, &(request[0]) ) );
550  CHK_MPI( MPI_Irecv( &(out_p_vec[recv_from_rank][0]),
551  (out_p_vec[recv_from_rank].size() * sizeof(T)),
552  MPI_BYTE,
553  recv_from_rank, 0, mpi_comm, &(request[1]) ) );
554  CHK_MPI( MPI_Waitall(2, request, MPI_STATUSES_IGNORE) );
555  }
556 }
MPI_Op mpi_typeof(std::logical_or< T >)
Definition: mpi.hpp:114
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the call graph for this function:

template<typename T , typename Partitioner >
void havoqgt::mpi::mpi_all_to_all_better ( std::vector< T > &  in_vec,
std::vector< T > &  out_vec,
Partitioner &  owner,
MPI_Comm  mpi_comm 
)

Definition at line 320 of file mpi.hpp.

321  {
322  int mpi_size(0), mpi_rank(0);
323  CHK_MPI( MPI_Comm_size( mpi_comm, &mpi_size) );
324  CHK_MPI( MPI_Comm_rank( mpi_comm, &mpi_rank) );
325 
326  std::vector<int> send_counts(mpi_size,0);
327  std::vector<int> send_disps(mpi_size,0);
328  std::vector<int> recv_counts(mpi_size,0);
329  std::vector<int> recv_disps(mpi_size,0);
330 
331  //sort send vector by owner
332  //std::sort(in_vec.begin(), in_vec.end(), owner_sort<T,Partitioner>(owner));
333 
334  //calc send counts
335  for(size_t i=0; i<in_vec.size(); ++i) {
336  send_counts[owner(in_vec[i], true)] += sizeof(T); //sizeof(t) lets us use PODS
337  //if(i>0) {
338  // assert(owner(inout_vec[i-1]) <= owner(inout_vec[i]));
339  //
340  }
341 
342 
343  //cacl send disps
344  std::partial_sum(send_counts.begin(), send_counts.end(), send_disps.begin());
345  for(size_t i=0; i<send_disps.size(); ++i) {
346  send_disps[i] -= send_counts[i]; //set to 0 offset
347  }
348 
349 
350  { //rearrange instead of sorting
351  std::vector<T> order_vec(in_vec.size());
352  std::vector<int> temp_arrange(mpi_size,0);
353  for(size_t i=0; i<in_vec.size(); ++i) {
354  int dest_rank = owner(in_vec[i], false);
355  assert (dest_rank >=0 && dest_rank < mpi_size);
356 
357  size_t dest_offset = send_disps[dest_rank] + temp_arrange[dest_rank];
358  temp_arrange[dest_rank] += sizeof(T);
359  dest_offset /= sizeof(T);
360  order_vec[dest_offset] = in_vec[i];
361  }
362  in_vec.swap(order_vec);
363  }
364 
365 
366  //exchange send counts
367  CHK_MPI( MPI_Alltoall( (void*) &(send_counts[0]), sizeof(int), MPI_BYTE,
368  (void*) &(recv_counts[0]), sizeof(int), MPI_BYTE,
369  mpi_comm ) );
370 
371 
372  //Allocate recv vector
373  int total_recv = std::accumulate(recv_counts.begin(), recv_counts.end(), 0);
374  out_vec.resize(total_recv/sizeof(T));
375 
376  //cacl recv disps
377  std::partial_sum(recv_counts.begin(), recv_counts.end(), recv_disps.begin());
378  for(size_t i=0; i<recv_disps.size(); ++i) {
379  recv_disps[i] -= recv_counts[i]; //set to 0 offset
380  }
381 
382  //perform actual alltoallv
383  void* send_ptr = in_vec.empty() ? NULL :(void*) &(in_vec[0]);
384  void* recv_ptr = out_vec.empty() ? NULL :(void*) &(out_vec[0]);
385  CHK_MPI( MPI_Alltoallv( send_ptr, &(send_counts[0]), &(send_disps[0]), MPI_BYTE,
386  recv_ptr, &(recv_counts[0]), &(recv_disps[0]), MPI_BYTE,
387  mpi_comm) );
388 }
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the caller graph for this function:

template<typename T >
void havoqgt::mpi::mpi_all_to_all_in_place ( std::vector< T > &  in_out_vec,
size_t  count,
MPI_Comm  mpi_comm 
)

Definition at line 397 of file mpi.hpp.

398  {
399  const int HAVOQGT_TAG=42;
400  int mpi_size(0), mpi_rank(0);
401  CHK_MPI( MPI_Comm_size( mpi_comm, &mpi_size) );
402  CHK_MPI( MPI_Comm_rank( mpi_comm, &mpi_rank) );
403  assert( in_out_vec.size() == count * mpi_size );
404 
405  for(int i=0; i<mpi_size; ++i) {
406  //start inner loop at i to avoid re-exchanging data
407  for(int j=i; j<mpi_size; ++j) {
408  if( mpi_rank == i) {
409  MPI_Status status;
410  CHK_MPI( MPI_Sendrecv_replace(&(in_out_vec[j*count]), count,
411  mpi_typeof(T()),
412  j, HAVOQGT_TAG,
413  j, HAVOQGT_TAG,
414  mpi_comm, &status) );
415 
416  } else if(mpi_rank == j) {
417  MPI_Status status;
418  CHK_MPI( MPI_Sendrecv_replace(&(in_out_vec[i*count]), count,
419  mpi_typeof(T()),
420  i, HAVOQGT_TAG,
421  i, HAVOQGT_TAG,
422  mpi_comm, &status) );
423 
424 
425  }
426  }
427  }
428 }
MPI_Op mpi_typeof(std::logical_or< T >)
Definition: mpi.hpp:114
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the call graph for this function:

template<typename T >
void havoqgt::mpi::mpi_bcast ( T &  data,
int  root,
MPI_Comm  comm 
)

Definition at line 559 of file mpi.hpp.

560 {
561  CHK_MPI( MPI_Bcast(&data, sizeof(data), MPI_BYTE, root, comm));
562 }
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the caller graph for this function:

int havoqgt::mpi::mpi_comm_rank ( )
inline

Definition at line 587 of file mpi.hpp.

587  {
588  int rank;
589  CHK_MPI( MPI_Comm_rank(MPI_COMM_WORLD, &rank) );
590  return rank;
591 }
#define CHK_MPI(a)
Definition: mpi.hpp:68
int havoqgt::mpi::mpi_comm_size ( )
inline

Definition at line 593 of file mpi.hpp.

593  {
594  int size;
595  CHK_MPI( MPI_Comm_size(MPI_COMM_WORLD, &size) );
596  return size;
597 }
#define CHK_MPI(a)
Definition: mpi.hpp:68
MPI_Datatype havoqgt::mpi::mpi_typeof ( char  )
inline

Definition at line 80 of file mpi.hpp.

80 {return MPI_CHAR;}

Here is the caller graph for this function:

MPI_Datatype havoqgt::mpi::mpi_typeof ( signed  short)
inline

Definition at line 81 of file mpi.hpp.

81 {return MPI_SHORT;}
MPI_Datatype havoqgt::mpi::mpi_typeof ( unsigned  char)
inline

Definition at line 84 of file mpi.hpp.

84 {return MPI_UNSIGNED_CHAR;}
MPI_Datatype havoqgt::mpi::mpi_typeof ( unsigned long long  )
inline

Definition at line 88 of file mpi.hpp.

88 {return MPI_UNSIGNED_LONG_LONG; }
MPI_Datatype havoqgt::mpi::mpi_typeof ( signed long  long)
inline

Definition at line 89 of file mpi.hpp.

89 {return MPI_LONG_LONG_INT;}
MPI_Datatype havoqgt::mpi::mpi_typeof ( double  )
inline

Definition at line 90 of file mpi.hpp.

90 {return MPI_DOUBLE;}
MPI_Datatype havoqgt::mpi::mpi_typeof ( long  double)
inline

Definition at line 91 of file mpi.hpp.

91 {return MPI_LONG_DOUBLE;}
MPI_Datatype havoqgt::mpi::mpi_typeof ( std::pair< int, int >  )
inline

Definition at line 92 of file mpi.hpp.

92 {return MPI_2INT;}
MPI_Datatype havoqgt::mpi::mpi_typeof ( std::pair< float, int >  )
inline

Definition at line 93 of file mpi.hpp.

93 {return MPI_FLOAT_INT;}
MPI_Datatype havoqgt::mpi::mpi_typeof ( std::pair< double, int >  )
inline

Definition at line 94 of file mpi.hpp.

94 {return MPI_DOUBLE_INT;}
MPI_Datatype havoqgt::mpi::mpi_typeof ( std::pair< long double, int >  )
inline

Definition at line 95 of file mpi.hpp.

95 {return MPI_LONG_DOUBLE_INT;}
MPI_Datatype havoqgt::mpi::mpi_typeof ( std::pair< short, int >  )
inline

Definition at line 96 of file mpi.hpp.

96 {return MPI_SHORT_INT;}
template<typename T >
MPI_Op havoqgt::mpi::mpi_typeof ( std::greater< T >  )

Definition at line 99 of file mpi.hpp.

99 { return MPI_MAX; }
template<typename T >
MPI_Op havoqgt::mpi::mpi_typeof ( std::less< T >  )

Definition at line 102 of file mpi.hpp.

102 { return MPI_MIN; }
template<typename T >
MPI_Op havoqgt::mpi::mpi_typeof ( std::plus< T >  )

Definition at line 105 of file mpi.hpp.

105 { return MPI_SUM; }
template<typename T >
MPI_Op havoqgt::mpi::mpi_typeof ( std::multiplies< T >  )

Definition at line 108 of file mpi.hpp.

108 { return MPI_PROD; }
template<typename T >
MPI_Op havoqgt::mpi::mpi_typeof ( std::logical_and< T >  )

Definition at line 111 of file mpi.hpp.

111 { return MPI_LAND; }
template<typename T >
MPI_Op havoqgt::mpi::mpi_typeof ( std::logical_or< T >  )

Definition at line 114 of file mpi.hpp.

114 { return MPI_LOR; }
void havoqgt::mpi::mpi_yield_barrier ( MPI_Comm  mpi_comm)

Definition at line 156 of file mpi.hpp.

156  {
157  MPI_Request request;
158  CHK_MPI(MPI_Ibarrier(mpi_comm, &request));
159 
160  for(; ;) {
161  MPI_Status status;
162  int is_done = false;
163  MPI_Test(&request, &is_done, &status);
164  if (is_done) {
165  return;
166  } else {
167  sched_yield();
168  }
169  }
170 }
#define CHK_MPI(a)
Definition: mpi.hpp:68

Here is the caller graph for this function:

template<typename SegmentManager >
bool havoqgt::mpi::operator!= ( const typename delegate_partitioned_graph< SegmentManager >::vertex_iterator &  x,
const typename delegate_partitioned_graph< SegmentManager >::vertex_iterator &  y 
)
inline

Definition at line 146 of file vertex_iterator.hpp.

147  {
148  return !(x.is_equal(y));
149 }

Here is the call graph for this function:

template<typename SegmentManager >
bool havoqgt::mpi::operator!= ( const typename delegate_partitioned_graph< SegmentManager >::edge_iterator &  x,
const typename delegate_partitioned_graph< SegmentManager >::edge_iterator &  y 
)
inline

Definition at line 149 of file edge_iterator.hpp.

150  {
151  return !(x.is_equal(y));
152 }

Here is the call graph for this function:

template<typename SegmentManager >
bool havoqgt::mpi::operator== ( const typename delegate_partitioned_graph< SegmentManager >::vertex_iterator &  x,
const typename delegate_partitioned_graph< SegmentManager >::vertex_iterator &  y 
)
inline

Definition at line 138 of file vertex_iterator.hpp.

139  {
140  return x.is_equal(y);
141 
142 }

Here is the call graph for this function:

template<typename SegmentManager >
bool havoqgt::mpi::operator== ( const typename delegate_partitioned_graph< SegmentManager >::edge_iterator &  x,
const typename delegate_partitioned_graph< SegmentManager >::edge_iterator &  y 
)
inline

Definition at line 141 of file edge_iterator.hpp.

142  {
143  return x.is_equal(y);
144 
145 }

Here is the call graph for this function:

template<typename TGraph , typename PRData >
void havoqgt::mpi::page_rank ( TGraph &  g,
PRData &  pr_data 
)

Definition at line 166 of file page_rank.hpp.

166  {
167  typedef pr_visitor<TGraph, PRData> visitor_type;
168  visitor_type::rank_data(&pr_data);
169  typedef visitor_queue< visitor_type, pr_queue, TGraph > visitor_queue_type;
170 
171  visitor_queue_type vq(&g);
172  vq.init_visitor_traversal();
173  pr_data.all_reduce();
174 }
template<typename TGraph , typename PathData , typename EdgeWeight >
void havoqgt::mpi::single_source_shortest_path ( TGraph &  g,
PathData &  path_data,
EdgeWeight &  edge_data,
typename TGraph::vertex_locator  s 
)

Definition at line 172 of file single_source_shortest_path.hpp.

175  {
176 
177  typedef sssp_visitor<TGraph, PathData, EdgeWeight> visitor_type;
178  visitor_type::set_path_data(&path_data);
179  visitor_type::set_edge_weight(&edge_data);
180  typedef visitor_queue< visitor_type, sssp_queue, TGraph > visitor_queue_type;
181 
182  visitor_queue_type vq(&g);
183  vq.init_visitor_traversal(s);
184 }
template<typename TGraph >
uint64_t havoqgt::mpi::triangle_count ( TGraph &  g,
typename TGraph::vertex_locator  s 
)

Definition at line 209 of file triangle_count.hpp.

209  {
210  typedef TGraph graph_type;
211  typedef typename mpi::triangle_count_visitor<TGraph> visitor_type;
212  visitor_type::m_ptr_graph = &g;
213  typedef visitor_queue< visitor_type, triangle_priority_queue, graph_type > visitor_queue_type;
214 
215  visitor_queue_type vq(&g);
216 
217  typename graph_type::template vertex_data<uint64_t, std::allocator<uint64_t> > tc_data(g);
218  visitor_type::set_tc_data(&tc_data);
219 
220  vq.init_visitor_traversal(s);
221 
222  return tc_data.global_accumulate();
223 }
hmpi::delegate_partitioned_graph< segment_manager_t > graph_type

Here is the caller graph for this function: