HavoqGT
gen_preferential_attachment_edge_list.hpp
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2013, Lawrence Livermore National Security, LLC.
3  * Produced at the Lawrence Livermore National Laboratory.
4  * Written by Roger Pearce <rpearce@llnl.gov>.
5  * LLNL-CODE-644630.
6  * All rights reserved.
7  *
8  * This file is part of HavoqGT, Version 0.1.
9  * For details, see https://computation.llnl.gov/casc/dcca-pub/dcca/Downloads.html
10  *
11  * Please also read this link – Our Notice and GNU Lesser General Public License.
12  * http://www.gnu.org/licenses/old-licenses/lgpl-2.1.html
13  *
14  * This program is free software; you can redistribute it and/or modify it under
15  * the terms of the GNU Lesser General Public License (as published by the Free
16  * Software Foundation) version 2.1 dated February 1999.
17  *
18  * This program is distributed in the hope that it will be useful, but WITHOUT ANY
19  * WARRANTY; without even the IMPLIED WARRANTY OF MERCHANTABILITY or FITNESS FOR A
20  * PARTICULAR PURPOSE. See the terms and conditions of the GNU General Public
21  * License for more details.
22  *
23  * You should have received a copy of the GNU Lesser General Public License along
24  * with this program; if not, write to the Free Software Foundation, Inc.,
25  * 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26  *
27  * OUR NOTICE AND TERMS AND CONDITIONS OF THE GNU GENERAL PUBLIC LICENSE
28  *
29  * Our Preamble Notice
30  *
31  * A. This notice is required to be provided under our contract with the
32  * U.S. Department of Energy (DOE). This work was produced at the Lawrence
33  * Livermore National Laboratory under Contract No. DE-AC52-07NA27344 with the DOE.
34  *
35  * B. Neither the United States Government nor Lawrence Livermore National
36  * Security, LLC nor any of their employees, makes any warranty, express or
37  * implied, or assumes any liability or responsibility for the accuracy,
38  * completeness, or usefulness of any information, apparatus, product, or process
39  * disclosed, or represents that its use would not infringe privately-owned rights.
40  *
41  * C. Also, reference herein to any specific commercial products, process, or
42  * services by trade name, trademark, manufacturer or otherwise does not
43  * necessarily constitute or imply its endorsement, recommendation, or favoring by
44  * the United States Government or Lawrence Livermore National Security, LLC. The
45  * views and opinions of authors expressed herein do not necessarily state or
46  * reflect those of the United States Government or Lawrence Livermore National
47  * Security, LLC, and shall not be used for advertising or product endorsement
48  * purposes.
49  *
50  */
51 
52 #ifndef HAVOQGT_MPI_GEN_PREFERENTIAL_ATTACHMENT_EDGE_LIST_HPP_INCLUDED
53 #define HAVOQGT_MPI_GEN_PREFERENTIAL_ATTACHMENT_EDGE_LIST_HPP_INCLUDED
54 
55 #include <havoqgt/mpi.hpp>
56 #include <havoqgt/detail/hash.hpp>
58 #include <boost/random.hpp>
59 #include <algorithm>
60 
61 namespace havoqgt { namespace mpi {
62 
63 
64 template <typename EdgeType>
65 void gen_preferential_attachment_edge_list(std::vector<EdgeType>& local_edges,
66  uint64_t in_base_seed,
67  uint64_t in_node_scale,
68  uint64_t in_edge_scale,
69  double in_beta,
70  double in_prob_rewire,
71  MPI_Comm in_comm)
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 }
238 
239 }} //end namespace havoqgt::mpi
240 
241 
242 #endif //end HAVOQGT_MPI_GEN_PREFERENTIAL_ATTACHMENT_EDGE_LIST_HPP_INCLUDED
243 
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
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)
#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