HavoqGT
generate_rmat.cpp
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 #include <boost/bind.hpp>
53 #include <boost/function.hpp>
58 #include <havoqgt/environment.hpp>
61 #include <iostream>
62 #include <assert.h>
63 #include <deque>
64 #include <utility>
65 #include <algorithm>
66 #include <functional>
67 #include <fstream> // std::ifstream
68 
69 
70 // notes for how to setup a good test
71 // take rank * 100 and make edges between (all local)
72 // Make one vert per rank a hub.
73 
74 using namespace havoqgt;
75 namespace hmpi = havoqgt::mpi;
76 using namespace havoqgt::mpi;
77 
78 void usage() {
79  if(havoqgt_env()->world_comm().rank() == 0) {
80  std::cerr << "Usage: -s <int> -d <int> -o <string>\n"
81  << " -s <int> - RMAT graph Scale (default 17)\n"
82  << " -d <int> - delegate threshold (Default is 1048576)\n"
83  << " -o <string> - output graph base filename\n"
84  << " -h - print help and exit\n\n";
85 
86  }
87 }
88 
89 void parse_cmd_line(int argc, char** argv, uint64_t& scale, uint64_t& delegate_threshold, std::string& output_filename) {
90  if(havoqgt_env()->world_comm().rank() == 0) {
91  std::cout << "CMD line:";
92  for (int i=0; i<argc; ++i) {
93  std::cout << " " << argv[i];
94  }
95  std::cout << std::endl;
96  }
97 
98  bool found_output_filename = false;
99  scale = 17;
100  delegate_threshold = 1048576;
101 
102  char c;
103  bool prn_help = false;
104  while ((c = getopt(argc, argv, "s:d:o:h ")) != -1) {
105  switch (c) {
106  case 'h':
107  prn_help = true;
108  break;
109  case 's':
110  scale = atoll(optarg);
111  break;
112  case 'd':
113  delegate_threshold = atoll(optarg);
114  break;
115  case 'o':
116  found_output_filename = true;
117  output_filename = optarg;
118  break;
119  default:
120  std::cerr << "Unrecognized option: "<<c<<", ignore."<<std::endl;
121  prn_help = true;
122  break;
123  }
124  }
125  if (prn_help || !found_output_filename) {
126  usage();
127  exit(-1);
128  }
129 }
130 
131 int main(int argc, char** argv) {
132 
134 
136 
137  int mpi_rank(0), mpi_size(0);
138 
139  havoqgt_init(&argc, &argv);
140  {
141  int mpi_rank = havoqgt_env()->world_comm().rank();
142  int mpi_size = havoqgt_env()->world_comm().size();
144 
145  if (mpi_rank == 0) {
146 
147  std::cout << "MPI initialized with " << mpi_size << " ranks." << std::endl;
149  }
151 
152  uint64_t num_vertices = 1;
153  uint64_t vert_scale;
154  uint64_t hub_threshold;
155  std::string fname_output;
156 
157  parse_cmd_line(argc, argv, vert_scale, hub_threshold, fname_output);
158 
159  num_vertices <<= vert_scale;
160  if (mpi_rank == 0) {
161  std::cout << "Building Graph500"<< std::endl
162  << "Building graph Scale: " << vert_scale << std::endl
163  << "Hub threshold = " << hub_threshold << std::endl
164  << "File name = " << fname_output << std::endl;
165  }
166 
167  havoqgt::distributed_db ddb(havoqgt::db_create(), fname_output.c_str());
168 
169  segment_manager_t* segment_manager = ddb.get_segment_manager();
170  bip::allocator<void, segment_manager_t> alloc_inst(segment_manager);
171 
172  //Generate RMAT graph
173  uint64_t num_edges_per_rank = num_vertices * 16 / mpi_size;
174  havoqgt::rmat_edge_generator rmat(uint64_t(5489) + uint64_t(mpi_rank) * 3ULL,
175  vert_scale, num_edges_per_rank,
176  0.57, 0.19, 0.19, 0.05, true, true);
177 
178 
179  if (mpi_rank == 0) {
180  std::cout << "Generating new graph." << std::endl;
181  }
182  graph_type *graph = segment_manager->construct<graph_type>
183  ("graph_obj")
184  (alloc_inst, MPI_COMM_WORLD, rmat, rmat.max_vertex_id(), hub_threshold);
185 
186 
188  if (mpi_rank == 0) {
189  std::cout << "Graph Ready, Calculating Stats. " << std::endl;
190  }
191 
192 
193 
194  for (int i = 0; i < mpi_size; i++) {
195  if (i == mpi_rank) {
196  double percent = double(segment_manager->get_free_memory()) /
197  double(segment_manager->get_size());
198  std::cout << "[" << mpi_rank << "] " << segment_manager->get_free_memory()
199  << "/" << segment_manager->get_size() << " = " << percent << std::endl;
200  }
202  }
203 
204  graph->print_graph_statistics();
205 
206 
207  //
208  // Calculate max degree
209  uint64_t max_degree(0);
210  for (auto citr = graph->controller_begin(); citr != graph->controller_end(); ++citr) {
211  max_degree = std::max(max_degree, graph->degree(*citr));
212  }
213 
214  uint64_t global_max_degree = havoqgt::mpi::mpi_all_reduce(max_degree, std::greater<uint64_t>(), MPI_COMM_WORLD);
215 
217 
218  if (mpi_rank == 0) {
219  std::cout << "Max Degree = " << global_max_degree << std::endl;
220  }
221 
223 
224  } //END Main MPI
226  return 0;
227 }
hmpi::delegate_partitioned_graph< segment_manager_t > graph_type
int main(int argc, char **argv)
void barrier() const
Definition: mpi.hpp:621
void havoqgt_finalize()
uint64_t degree(vertex_locator locator) const
Returns the degree of a vertex.
const communicator & world_comm() const
int size() const
Definition: mpi.hpp:618
T mpi_all_reduce(T in_d, Op in_op, MPI_Comm mpi_comm)
Definition: mpi.hpp:176
int rank() const
Definition: mpi.hpp:619
void usage()
environment * havoqgt_env()
mapped_type::segment_manager segment_manager_type
segment_manager_type * get_segment_manager()
old_environment & get_environment()
void parse_cmd_line(int argc, char **argv, uint64_t &scale, uint64_t &delegate_threshold, std::string &output_filename)
havoqgt::distributed_db::segment_manager_type segment_manager_t
void havoqgt_init(int *argc, char ***argv)