graph.hh

Go to the documentation of this file.
1 
6 #ifndef graph_hh
7 #define graph_hh
8 
9 #include "basics/exceptions.hh"
10 #include "basics/debug.hh"
11 #include "toolbox/sequence.hh"
12 #ifdef HAS_MPI
13 #include <mpi.h>
14 #endif
15 
16 #ifdef HAS_METIS
17 #include <metis.h>
18 #endif
19 
20 #define GRAPH_OUTPUT 0
21 
22 
23 namespace concepts
24 {
25 
26  // *********************************************************** GraphVertex **
27 
28 
31  template<class F>
33  {
34  public:
36  GraphVertex(const F Value, const int Weight=1, const int Domain=0) : Value_(Value), Weight_(Weight), Domain_(Domain) {}
38  virtual void setValue(const F Value);
39  virtual void setWeight(const int Weight);
40  virtual void addWeight(const int Weight);
41  virtual void setNeighbours(const Sequence<F> & Neighbours);
42  virtual void addNeighbour(const F Value);
43  virtual void addWeightedNeighbour(const F Value, const int Weight = 1);
44  virtual void setDomain(const int Domain);
45  virtual void clearNeighbours();
46  virtual F getValue() const;
47  virtual int getWeight() const;
48  virtual int getDomain() const;
49  virtual Sequence<F> getNeighbours() const;
50  virtual std::ostream& info(std::ostream& os) const;
51  protected:
52  F Value_;
54  int Weight_;
55  int Domain_;
56  };
57 
58  template<class F>
59  std::ostream& GraphVertex<F>::info(std::ostream& os) const
60  {
61  return os << concepts::typeOf(*this)<<": (" << this->Value_
62  << "," << this->Neighbours_
63  << "," << this->Weight_
64  << "," << this->Domain_ << ")" ;
65  }
66 
67  template<class F>
68  std::ostream& operator<<(std::ostream& os, const GraphVertex<F> & Value)
69  {
70  return Value.info(os);
71  }
72 
73  // ***************************************************************** Graph **
74 
75 
78  template<class F>
79  class Graph
80  {
81  public:
82  Graph() {}
83  ~Graph() {Vertexes_.clear();}
84  virtual void addVertex(const F Value, const int Weight = 1);
85  virtual void addEdge(const F Value1, const F Value2);
86  virtual void addWeightedEdge(const F Value1, const F Value2,
87  const int Weight);
89  virtual void GetCompressedFormat(Sequence<int> & Cumuled, Sequence<int> & CNeighbours, Sequence<int> & Weights) const;
90 #ifdef HAS_METIS
91  virtual Sequence<int> GetMetisSubdomains(const int NbParts = 2) const;
92 #endif
93  virtual void SetSubdomains(const Sequence<int> & Subdomains);
94  virtual void SetSubdomains(const int NbParts = 2);
95  void EmptyGraph() {Vertexes_.clear(); };
96  virtual Sequence<F> GetSubdomain(const int Part) const;
97  virtual Sequence<bool> GetSubdomainBool(const int Part) const;
98  virtual std::ostream& info(std::ostream& os) const;
99  protected:
101  };
102 
103  template<class F>
104  std::ostream& Graph<F>::info(std::ostream& os) const
105  {
106  int n=1;
107  for(typename std::vector<GraphVertex<F> >::const_iterator i
108  = this->Vertexes_.begin(); i != this->Vertexes_.end(); ++i, ++n)
109  os << "Point " << n << "(value " << i->getValue() << ") belongs to graph " << i->getDomain() << std::endl;
110 
111  return os;
112  }
113 
114  template<class F>
115  std::ostream& operator<<(std::ostream& os, const Graph<F> & Value)
116  {
117  return Value.info(os);
118  }
119 
120  template<class F>
121  void GraphVertex<F>::setValue(const F Value)
122  {
123  this->Value_ = Value;
124  }
125 
126  template<class F>
127  void GraphVertex<F>::setWeight(const int Weight)
128  {
129  this->Weight_ = Weight;
130  }
131 
132  template<class F>
133  void GraphVertex<F>::addWeight(const int Weight)
134  {
135  this->Weight_ += Weight;
136  }
137 
138  template<class F>
140  {
141  this->Neighbours_ = Neighbours;
142  }
143 
144  template<class F>
145  void GraphVertex<F>::addNeighbour(const F Value)
146  {
147  if ( this->Neighbours_.exist(Value) )
148  return;
149  this->Neighbours_.push_back(Value);
150  }
151 
152  template<class F>
153  void GraphVertex<F>::addWeightedNeighbour(const F Value, const int Weight)
154  {
155  if ( this->Neighbours_.exist(Value) )
156  return;
157  this->Neighbours_.push_back(Value);
158  this->Weight_ += Weight;
159  }
160 
161  template<class F>
162  void GraphVertex<F>::setDomain(const int Domain)
163  {
164  this->Domain_ = Domain;
165  }
166 
167  template<class F>
169  {
170  this->Neighbours_.clear();
171  }
172 
173  template<class F>
175  {
176  return this->Value_;
177  }
178 
179  template<class F>
181  {
182  return this->Weight_;
183  }
184 
185  template<class F>
187  {
188  return this->Neighbours_;
189  }
190 
191  template<class F>
193  {
194  return this->Domain_;
195  }
196 
197 
198 
203  template<class F>
204  void Graph<F>::addVertex(const F Value, const int Weight)
205  {
206 #ifdef DEBUG
207  for(typename std::vector<GraphVertex<F> >::const_iterator i
208  = this->Vertexes_.begin(); i != this->Vertexes_.end(); ++i)
209  if (i->getValue() == Value)
210  return;
211 #endif
212  GraphVertex<F> VValue(Value, Weight);
213  this->Vertexes_.push_back(VValue);
214  }
215 
224  template<class F>
225  void Graph<F>::addEdge(const F Value1, const F Value2)
226  {
227  bool Value1NotDone = true, Value2NotDone = true;
228  for(typename std::vector<GraphVertex<F> >::iterator i
229  = this->Vertexes_.begin();
230  i != this->Vertexes_.end() && (Value1NotDone || Value2NotDone); ++i)
231  {
232  if (i->getValue() == Value1)
233  {
234  Value1NotDone = false;
235  i->addNeighbour(Value2);
236  }
237  if (i->getValue() == Value2)
238  {
239  Value2NotDone = false;
240  i->addNeighbour(Value1);
241  }
242  }
243  }
244 
254  template<class F>
255  void Graph<F>::addWeightedEdge(const F Value1, const F Value2,
256  const int Weight)
257  {
258  bool Value1NotDone = true, Value2NotDone = true;
259  for(typename std::vector<GraphVertex<F> >::iterator i
260  = this->Vertexes_.begin();
261  i != this->Vertexes_.end() && (Value1NotDone || Value2NotDone); ++i)
262  {
263  if (i->getValue() == Value1)
264  {
265  Value1NotDone = false;
266  i->addWeightedNeighbour(Value2, Weight);
267  }
268  if (i->getValue() == Value2)
269  {
270  Value2NotDone = false;
271  i->addWeightedNeighbour(Value1, Weight);
272  }
273  }
274  }
275 
276 
277  template<class F>
279  {
280  Sequence<F> Result;
281  Result.clear();
282  for(typename std::vector<GraphVertex<F> >::const_iterator i
283  = this->Vertexes_.begin(); i != this->Vertexes_.end(); ++i)
284  {
285  Result = Result + i->getNeighbours();
286  }
287  return Result;
288  }
289 
290  template<class F>
291  void Graph<F>::GetCompressedFormat(Sequence<int> & Cumuled, Sequence<int> & CNeighbours, Sequence<int> & Weights) const
292  {
293  Cumuled.clear();
294  CNeighbours.clear();
295  Weights.clear();
296  Cumuled.push_back(0); // Cumuled(0) = 0;
297  Sequence<F> TempCNeighbours;
298  TempCNeighbours.clear();
299  int cml = 0;
300  int n = 0;
301  for(typename std::vector<GraphVertex<F> >::const_iterator i
302  = this->Vertexes_.begin(); i != this->Vertexes_.end(); ++i, ++n)
303  {
304  TempCNeighbours = TempCNeighbours + i->getNeighbours();
305  cml += i->getNeighbours().size();
306  Cumuled.push_back(cml);
307  //Cumuled(n+1) = Cumuled(n) + i->getNeighbours().size();
308  Weights.push_back(i->getWeight()); // Weights(n) = i->getWeight();
309 
310  }
311 
312 
313  for(typename std::vector<F>::const_iterator i
314  = TempCNeighbours.begin(); i != TempCNeighbours.end(); ++i)
315  {
316  /* We retrieve compressed position of the neighbour (i.e. the position in the vertex
317  sequence, starting from n=0 for the first position */
318  n = 0;
319  typename std::vector<GraphVertex<F> >::const_iterator j = this->Vertexes_.begin();
320  for(; (j != this->Vertexes_.end()) && (*i != j->getValue()); ++j,
321  ++n);
322  CNeighbours.push_back(n);
323  }
324 
325  TempCNeighbours.clear();
326 
327  }
328 
329 #ifdef HAS_METIS
330  template<class F>
331  Sequence<int> Graph<F>::GetMetisSubdomains(const int NbParts) const
332  {
333  Sequence<int> Subdomains;
334  Subdomains.clear();
335  Subdomains.resize(this->Vertexes_.size());
336  concepts::Sequence<int> Cumuled;
337  concepts::Sequence<int> CNeighbours;
338  concepts::Sequence<int> Weights;
339  this->GetCompressedFormat(Cumuled, CNeighbours, Weights);
340 
341  idx_t nvtxs = idx_t(Cumuled.size()-1);
342  idx_t ncon = 1;
343  idx_t *xadj = new idx_t[Cumuled.size()];
344 
345  int n=0;
346  for(typename std::vector<int>::const_iterator i
347  = Cumuled.begin(); i != Cumuled.end(); ++i, ++n)
348  {
349  xadj[n] = *i;
350  }
351 
352  idx_t * vsize = NULL;
353 
354  idx_t * adjncy = new idx_t[CNeighbours.size()];
355  n=0;
356  for(typename std::vector<int>::const_iterator i
357  = CNeighbours.begin(); i != CNeighbours.end(); ++i, ++n)
358  {
359  adjncy[n] = *i;
360  }
361 
362  idx_t * vwgt = new idx_t[Weights.size()];
363  n=0;
364  for(typename std::vector<int>::const_iterator i
365  = Weights.begin(); i != Weights.end(); ++i, ++n)
366  {
367  vwgt[n] = *i;
368  }
369 
370  idx_t *adjwgt = NULL;
371  idx_t nparts = idx_t(NbParts);
372  real_t *tpwgts = NULL, *ubvec= NULL;
373  idx_t *options=NULL;
374  idx_t objval;
375  int vertex_size = this->Vertexes_.size();
376  idx_t* part = new idx_t[vertex_size];
377 
378  int result = METIS_PartGraphKway( &nvtxs, &ncon, xadj, adjncy, vwgt,
379  vsize, adjwgt, &nparts, tpwgts,
380  ubvec, options, &objval, part);
381 
382  conceptsAssert3(result==METIS_OK, Assertion(), "Partition graph failed");
383 #ifdef HAS_MPI
384  int mpisize; MPI_Comm_size(MPI_COMM_WORLD, &mpisize);
385  int mpirank; MPI_Comm_rank(MPI_COMM_WORLD, &mpirank);
386  DEBUGL(GRAPH_OUTPUT, "Converting vector of size " << vertex_size << " on process " << mpirank);
387  int* mpipart = new int[vertex_size];
388  for(n=0 ; n < vertex_size ; n++)
389  {
390  mpipart[n] = int(part[n]);
391  }
392  DEBUGL(GRAPH_OUTPUT, "Starting Broadcasting on process " << mpirank);
393  if (mpisize > 1)
394  {
395  MPI_Barrier(MPI_COMM_WORLD);
396  MPI_Bcast(mpipart, vertex_size, MPI_INT, 0, MPI_COMM_WORLD);
397  MPI_Barrier(MPI_COMM_WORLD);
398  }
399  DEBUGL(GRAPH_OUTPUT, "Broadcast ended on process " << mpirank);
400 #endif
401 
402  n=0;
403  for(typename std::vector<int>::iterator i
404  = Subdomains.begin(); i!= Subdomains.end(); ++i, ++n)
405  {
406 #ifdef HAS_MPI
407  *i = mpipart[n];
408 #else
409  *i = part[n];
410 #endif
411  }
412 
413  delete[] xadj;
414  delete[] adjncy;
415  delete[] vwgt;
416  delete[] part;
417 #ifdef HAS_MPI
418  delete[] mpipart;
419 #endif
420  Cumuled.clear();
421  CNeighbours.clear();
422  Weights.clear();
423 
424  return Subdomains;
425  }
426 
427 #endif
428 
429  template<class F>
430  void Graph<F>::SetSubdomains(const Sequence<int> & Subdomains)
431  {
432 
433  assert(Subdomains.size() == this->Vertexes_.size());
434  typename std::vector<GraphVertex<F> >::iterator iv =
435  this->Vertexes_.begin();
436  typename std::vector<int>::const_iterator i = Subdomains.begin();
437  for( ; i!=Subdomains.end(); ++iv, ++i)
438  {
439  iv->setDomain(*i);
440  }
441  }
442 
443 
444  template<class F>
445  void Graph<F>::SetSubdomains(const int NbParts)
446  {
447  if (NbParts == 1)
448  {
449  Sequence<int> Subdomains;
450  Subdomains.clear();
451  Subdomains.resize(this->Vertexes_.size());
452  for(typename std::vector<int>::iterator i
453  = Subdomains.begin(); i!= Subdomains.end(); ++i)
454  *i = 0;
455  this->SetSubdomains(Subdomains);
456  Subdomains.clear();
457  }
458  else
459  {
460 #ifdef HAS_METIS
461  this->SetSubdomains(this->GetMetisSubdomains(NbParts));
462 #else
463  // Give a dummy graph splitting implementation
464  // This implementation is not optimized, but works
465  Sequence<int> Subdomains;
466  Subdomains.clear();
467  Subdomains.resize(this->Vertexes_.size());
468  int n=0;
469  for(typename std::vector<int>::iterator i
470  = Subdomains.begin(); i!= Subdomains.end(); ++i, ++n)
471  *i = (n * NbParts) / this->Vertexes_.size();
472  this->SetSubdomains(Subdomains);
473  Subdomains.clear();
474 #endif
475  }
476  }
477 
478  template<class F>
479  Sequence<F> Graph<F>::GetSubdomain(const int Part) const
480  {
481  Sequence<F> ListVertexes;
482  ListVertexes.clear();
483  typename std::vector<GraphVertex<F> >::const_iterator iv =
484  this->Vertexes_.begin();
485  for ( ; iv != this->Vertexes_.end(); ++iv)
486  {
487  if ( iv->getDomain() == Part )
488  {
489  ListVertexes.push_back(iv->getValue());
490  }
491  }
492  return ListVertexes;
493  }
494 
495  template<class F>
497  {
498  Sequence<bool> ListVertexes;
499  ListVertexes.clear();
500  typename std::vector<GraphVertex<F> >::const_iterator iv =
501  this->Vertexes_.begin();
502  for ( ; iv != this->Vertexes_.end(); ++iv)
503  {
504  if ( iv->getDomain() == Part )
505  ListVertexes.push_back(true);
506  else
507  ListVertexes.push_back(false);
508  }
509  return ListVertexes;
510  }
511 
512 }
513 
514 #endif // graph_hh
virtual std::ostream & info(std::ostream &os) const
Definition: graph.hh:59
virtual Sequence< bool > GetSubdomainBool(const int Part) const
Definition: graph.hh:496
virtual void GetCompressedFormat(Sequence< int > &Cumuled, Sequence< int > &CNeighbours, Sequence< int > &Weights) const
Definition: graph.hh:291
void EmptyGraph()
Definition: graph.hh:95
virtual void addWeightedEdge(const F Value1, const F Value2, const int Weight)
Add edge into graph (with weight modification)
Definition: graph.hh:255
virtual void clearNeighbours()
Definition: graph.hh:168
virtual void addVertex(const F Value, const int Weight=1)
Add vertex into graph.
Definition: graph.hh:204
virtual void setWeight(const int Weight)
Definition: graph.hh:127
virtual void setDomain(const int Domain)
Definition: graph.hh:162
virtual void addEdge(const F Value1, const F Value2)
Add edge into graph.
Definition: graph.hh:225
virtual void setValue(const F Value)
Definition: graph.hh:121
#define DEBUGL(doit, msg)
Template class to define a graph vertex.
Definition: graph.hh:33
virtual void SetSubdomains(const int NbParts=2)
Definition: graph.hh:445
Sequence< GraphVertex< F > > Vertexes_
Definition: graph.hh:100
virtual int getDomain() const
Definition: graph.hh:192
GraphVertex(const F Value, const int Weight=1, const int Domain=0)
Definition: graph.hh:36
virtual void addWeight(const int Weight)
Definition: graph.hh:133
virtual Sequence< F > getCompressedNeighbours() const
Definition: graph.hh:278
std::ostream & operator<<(std::ostream &os, const Level< dim > &c)
virtual Sequence< F > getNeighbours() const
Definition: graph.hh:186
Sequence with operations, output operator, and method of the particular element types.
Definition: sequence.hh:39
virtual F getValue() const
Definition: graph.hh:174
virtual void addNeighbour(const F Value)
Definition: graph.hh:145
virtual int getWeight() const
Definition: graph.hh:180
virtual std::ostream & info(std::ostream &os) const
Definition: graph.hh:104
virtual void setNeighbours(const Sequence< F > &Neighbours)
Definition: graph.hh:139
Sequence< F > Neighbours_
Definition: graph.hh:53
#define conceptsAssert3(cond, exc, msg)
Assert that a certain condition is fulfilled.
Definition: exceptions.hh:442
std::string typeOf(const T &t)
Return the typeid name of a class object.
Definition: output.hh:43
virtual Sequence< F > GetSubdomain(const int Part) const
Definition: graph.hh:479
#define GRAPH_OUTPUT
Definition: graph.hh:20
virtual void addWeightedNeighbour(const F Value, const int Weight=1)
Definition: graph.hh:153
virtual void SetSubdomains(const Sequence< int > &Subdomains)
Definition: graph.hh:430
Basic namespace for Concepts-2.
Definition: pml_formula.h:16
Template class to define a graph.
Definition: graph.hh:80
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich