HavoqGT
delegate_partitioned_graph.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_DELEGATE_PARTITIONED_GRAPH_HPP_INCLUDED
53 #define HAVOQGT_MPI_DELEGATE_PARTITIONED_GRAPH_HPP_INCLUDED
54 
55 
56 #include <limits>
57 #include <utility>
58 #include <stdint.h>
59 #include <functional>
60 
61 #include <boost/unordered_set.hpp>
62 #include <boost/unordered_map.hpp>
63 #include <boost/container/map.hpp>
64 #include <boost/interprocess/containers/map.hpp>
65 #include <boost/interprocess/containers/vector.hpp>
66 #include <boost/interprocess/allocators/allocator.hpp>
67 
68 #include <havoqgt/mpi.hpp>
69 #include <havoqgt/utilities.hpp>
74 
75 #ifdef DEBUG_DPG
76  #warning Debug MACRO is for delegate_partitioned_graph.
77  #define IS_DEBUGING true
78 #else
79  #define IS_DEBUGING false
80 #endif
81 
82 namespace havoqgt {
83 namespace mpi {
84 
85 namespace bip = boost::interprocess;
86 
97 template <typename SegementManager>
99  public:
100  template<typename T>
101  using SegmentAllocator = bip::allocator<T, SegementManager>;
102 
104  // Class Objects
107  class vertex_locator;
109  class edge_iterator;
111  class vertex_iterator;
113  template <typename T, typename Allocator>
114  class vertex_data;
116  template <typename T, typename SegManagerOther>
117  class edge_data;
119  class vert_info;
120 
123 
125  // Public Member Functions
128  template <typename Container>
129  delegate_partitioned_graph(const SegmentAllocator<void>& seg_allocator,
130  MPI_Comm mpi_comm,
131  Container& edges, uint64_t max_vertex,
132  uint64_t delegate_degree_threshold,
133  ConstructionState stop_after = GraphReady);
134 
135  template <typename Container>
136  void complete_construction(const SegmentAllocator<void>& seg_allocator,
137  MPI_Comm mpi_comm, Container& edges);
138  void print_graph_statistics();
139 
141  uint64_t locator_to_label(vertex_locator locator) const;
142 
144  vertex_locator label_to_locator(uint64_t label) const;
145 
147  edge_iterator edges_begin(vertex_locator locator) const;
148 
150  edge_iterator edges_end(vertex_locator locator) const;
151 
153  uint64_t degree(vertex_locator locator) const;
154 
156  uint64_t local_degree(vertex_locator locator) const;
157 
159  vertex_iterator vertices_begin() const;
160 
162  vertex_iterator vertices_end() const;
163 
165  bool is_label_delegate(uint64_t label) const;
166 
168  template <typename T, typename SegManagerOther>
169  vertex_data<T, SegManagerOther>* create_vertex_data(
170  SegManagerOther*,
171  const char *obj_name = nullptr) const;
172 
174  template <typename T, typename SegManagerOther>
175  vertex_data<T, SegManagerOther>* create_vertex_data(
176  const T& init, SegManagerOther*,
177  const char *obj_name = nullptr) const;
178 
180  template <typename T, typename SegManagerOther>
181  edge_data<T, SegManagerOther>* create_edge_data(
182  SegManagerOther*,
183  const char *obj_name = nullptr) const;
184 
186  template <typename T, typename SegManagerOther>
187  edge_data<T, SegManagerOther>* create_edge_data(
188  const T& init, SegManagerOther*,
189  const char *obj_name = nullptr) const;
190 
191  size_t num_local_vertices() const {
192  return m_owned_info.size();
193  }
194 
195  uint64_t max_global_vertex_id() {
196  return m_global_max_vertex;
197  }
198 
199  uint64_t max_local_vertex_id() {
200  return m_max_vertex;
201  }
202 
203 
204  size_t num_delegates() const {
205  return m_delegate_degree.size();
206  }
207 
208  uint32_t master(const vertex_locator& locator) const {
209  return locator.m_local_id % m_mpi_size;
210  }
211 
212  typedef typename bip::vector<vertex_locator, SegmentAllocator<vertex_locator> >
213  ::const_iterator controller_iterator;
214 
215  controller_iterator controller_begin() const {
216  return m_controller_locators.begin();
217  }
218 
219  controller_iterator controller_end() const {
220  return m_controller_locators.end();
221  }
222 
224  return *this==*b;
225  }
226 
228  if (m_mpi_size != other.m_mpi_size)
229  return false;
230  if (m_mpi_rank != other.m_mpi_rank)
231  return false;
232  if(m_owned_info != other.m_owned_info)
233  return false;
234  // if(m_owned_targets != other.m_owned_targets)
235  // return false;
236  if(m_delegate_info != other.m_delegate_info)
237  return false;
239  return false;
241  return false;
242  // if(m_delegate_targets != other.m_delegate_targets)
243  // return false;
245  return false;
247  return false;
249  return false;
251  return false;
253  return false;
254 
255  return true;
256  }
257 
259  other) {
260  return !(*this == other);
261  }
262 
263  private:
265  // Private Member Functions
268  void sync_global_hub_set(const boost::unordered_set<uint64_t>& local_hubs,
269  boost::unordered_set<uint64_t>& global_hubs,
270  bool local_change);
271 
272 
273  void initialize_low_meta_data(boost::unordered_set<uint64_t>& global_hub_set);
274  void initialize_high_meta_data(boost::unordered_set<uint64_t>& global_hubs);
275 
276  void initialize_edge_storage(const SegmentAllocator<void>& seg_allocator);
277 
278  template <typename Container>
279  void partition_low_degree(Container& unsorted_edges);
280 
281  template <typename InputIterator>
282  void count_high_degree_edges(InputIterator unsorted_itr,
283  InputIterator unsorted_itr_end,
284  boost::unordered_set<uint64_t>& global_hub_set);
285 
286 
287  template <typename Container>
288  void partition_high_degree(Container& unsorted_edges,
289  std::map< uint64_t, std::deque<OverflowSendInfo> > &transfer_info);
290 
291 
292  template <typename InputIterator>
293  void count_edge_degrees(InputIterator unsorted_itr,
294  InputIterator unsorted_itr_end,
295  boost::unordered_set<uint64_t>& global_hub_set,
296  uint64_t delegate_degree_threshold);
297 
298 
299  void send_high_info(std::vector< boost::container::map< uint64_t, uint64_t> >&
300  maps_to_send, int maps_to_send_element_count);
301 
302 
303  void send_vertex_info(uint64_t &high_vertex_count,
304  uint64_t delegate_degree_threshold,
305  std::vector< boost::container::map<
306  int, std::pair<uint64_t, uint64_t> > >& maps_to_send,
307  int maps_to_send_element_count);
308 
309 
310  void calculate_overflow(std::map< uint64_t, std::deque<OverflowSendInfo> >
311  &transfer_info);
312 
313  void generate_send_list(std::vector<uint64_t> &send_list, uint64_t num_send,
314  int send_id,
315  std::map< uint64_t, std::deque<OverflowSendInfo> > &transfer_info);
316 
317  void flush_graph();
318 
320  // Protected Data Members
322  protected:
323  uint64_t edge_chunk_size;
328  MPI_Comm m_mpi_comm;
329 
331 
332 
333  uint64_t m_max_vertex {0};
334  uint64_t m_global_max_vertex {0};
335  uint64_t m_global_edge_count {0};
336  uint64_t m_edges_high_count {0};
337  uint64_t m_edges_low_count {0};
338 
339  bip::vector<uint32_t, SegmentAllocator<uint32_t> > m_local_outgoing_count;
340  bip::vector<uint32_t, SegmentAllocator<uint32_t> > m_local_incoming_count;
341 
342  bip::vector<vert_info, SegmentAllocator<vert_info>> m_owned_info;
343  bip::vector<uint32_t, SegmentAllocator<uint32_t>> m_owned_info_tracker;
344  //bip::vector<vertex_locator, SegmentAllocator<vertex_locator>> m_owned_targets;
345  bip::offset_ptr<vertex_locator> m_owned_targets;
347 
348  // Delegate Storage
350 
351  bip::vector< uint64_t, SegmentAllocator<uint64_t> > m_delegate_info;
352  bip::vector< uint64_t, SegmentAllocator<uint64_t> > m_delegate_degree;
353  bip::vector< uint64_t, SegmentAllocator<uint64_t> > m_delegate_label;
354  // bip::vector< vertex_locator, SegmentAllocator<vertex_locator> >
355  // m_delegate_targets;
356  bip::offset_ptr<vertex_locator> m_delegate_targets;
358 
359  //Note: BIP only contains a map, not an unordered_map object.
360  /*boost::interprocess::unordered_map<
361  uint64_t, vertex_locator, boost::hash<uint64_t>, std::equal_to<uint64_t>,
362  SegmentAllocator< std::pair<uint64_t,vertex_locator> >
363  > m_map_delegate_locator;
364  */
365  bip::map<uint64_t, vertex_locator, std::less<uint64_t>, SegmentAllocator< std::pair<const uint64_t,vertex_locator> > > m_map_delegate_locator;
366 
367  bip::vector<vertex_locator, SegmentAllocator<vertex_locator> >
369 }; // class delegate_partitioned_graph
370 
371 
372 
373 } // namespace mpi
374 } // namespace havoqgt
375 
376 
377 #include <havoqgt/impl/log_step.hpp>
384 
386 
387 #endif //HAVOQGT_MPI_DELEGATE_PARTITIONED_GRAPH_HPP_INCLUDED
bool operator!=(delegate_partitioned_graph< SegementManager > &other)
bip::vector< vertex_locator, SegmentAllocator< vertex_locator > > m_controller_locators
edge_iterator edges_end(vertex_locator locator) const
Returns an end iterator for edges of a vertex.
void count_edge_degrees(InputIterator unsorted_itr, InputIterator unsorted_itr_end, boost::unordered_set< uint64_t > &global_hub_set, uint64_t delegate_degree_threshold)
bip::vector< vertex_locator, SegmentAllocator< vertex_locator > >::const_iterator controller_iterator
vertex_iterator vertices_end() const
Returns an end iterator for all local vertices.
void partition_low_degree(Container &unsorted_edges)
bip::vector< uint64_t, SegmentAllocator< uint64_t > > m_delegate_info
bip::vector< uint32_t, SegmentAllocator< uint32_t > > m_local_outgoing_count
uint64_t degree(vertex_locator locator) const
Returns the degree of a vertex.
void send_high_info(std::vector< boost::container::map< uint64_t, uint64_t > > &maps_to_send, int maps_to_send_element_count)
bip::vector< uint64_t, SegmentAllocator< uint64_t > > m_delegate_label
uint64_t locator_to_label(vertex_locator locator) const
Converts a vertex_locator to the vertex label.
edge_iterator edges_begin(vertex_locator locator) const
Returns a begin iterator for edges of a vertex.
bip::map< uint64_t, vertex_locator, std::less< uint64_t >, SegmentAllocator< std::pair< const uint64_t, vertex_locator > > > m_map_delegate_locator
bool is_label_delegate(uint64_t label) const
Tests if vertex label is a delegate.
uint64_t local_degree(vertex_locator locator) const
Returns the local degree of a vertex.
vertex_iterator vertices_begin() const
Returns a begin iterator for all local vertices.
void send_vertex_info(uint64_t &high_vertex_count, uint64_t delegate_degree_threshold, std::vector< boost::container::map< int, std::pair< uint64_t, uint64_t > > > &maps_to_send, int maps_to_send_element_count)
void complete_construction(const SegmentAllocator< void > &seg_allocator, MPI_Comm mpi_comm, Container &edges)
Edge storage allocation phase of graph construction.
uint32_t master(const vertex_locator &locator) const
void initialize_high_meta_data(boost::unordered_set< uint64_t > &global_hubs)
void partition_high_degree(Container &unsorted_edges, std::map< uint64_t, std::deque< OverflowSendInfo > > &transfer_info)
bool compare(delegate_partitioned_graph< SegementManager > *b)
bip::vector< vert_info, SegmentAllocator< vert_info > > m_owned_info
bip::vector< uint32_t, SegmentAllocator< uint32_t > > m_local_incoming_count
void generate_send_list(std::vector< uint64_t > &send_list, uint64_t num_send, int send_id, std::map< uint64_t, std::deque< OverflowSendInfo > > &transfer_info)
void sync_global_hub_set(const boost::unordered_set< uint64_t > &local_hubs, boost::unordered_set< uint64_t > &global_hubs, bool local_change)
Synchronizes hub set amongst all processes.
void initialize_edge_storage(const SegmentAllocator< void > &seg_allocator)
void calculate_overflow(std::map< uint64_t, std::deque< OverflowSendInfo > > &transfer_info)
vertex_locator label_to_locator(uint64_t label) const
Converts a vertex label to a vertex_locator.
bool operator==(delegate_partitioned_graph< SegementManager > &other)
bip::offset_ptr< vertex_locator > m_delegate_targets
void initialize_low_meta_data(boost::unordered_set< uint64_t > &global_hub_set)
bip::vector< uint64_t, SegmentAllocator< uint64_t > > m_delegate_degree
bip::vector< uint32_t, SegmentAllocator< uint32_t > > m_owned_info_tracker
bip::allocator< T, SegementManager > SegmentAllocator
delegate_partitioned_graph(const SegmentAllocator< void > &seg_allocator, MPI_Comm mpi_comm, Container &edges, uint64_t max_vertex, uint64_t delegate_degree_threshold, ConstructionState stop_after=GraphReady)
Constructor that initializes given and unsorted sequence of edges.
void count_high_degree_edges(InputIterator unsorted_itr, InputIterator unsorted_itr_end, boost::unordered_set< uint64_t > &global_hub_set)
bip::offset_ptr< vertex_locator > m_owned_targets
vertex_data< T, SegManagerOther > * create_vertex_data(SegManagerOther *, const char *obj_name=nullptr) const
Creates vertex_data of type T.
edge_data< T, SegManagerOther > * create_edge_data(SegManagerOther *, const char *obj_name=nullptr) const
Creates edge_data of type T.