OpenFPM_pdata  1.1.0
Project that contain the implementation of distributed structures
 All Data Structures Namespaces Functions Variables Typedefs Enumerations Friends Pages
metis_util.hpp
1 /*
2  * metis_util.hpp
3  *
4  * Created on: Nov 21, 2014
5  * Author: Pietro Incardona
6  */
7 
8 #ifndef METIS_UTIL_HPP
9 #define METIS_UTIL_HPP
10 
11 #include <iostream>
12 #include "metis.h"
13 #include "SubdomainGraphNodes.hpp"
14 #include "VTKWriter/VTKWriter.hpp"
15 
22 {
24  idx_t * nvtxs;
25 
30  idx_t * ncon;
31 
33  idx_t * xadj;
34 
36  idx_t * adjncy;
37 
39  idx_t * vwgt;
40 
42  idx_t * vsize;
44  idx_t * adjwgt;
46  idx_t * nparts;
48  real_t * tpwgts;
50  real_t * ubvec;
52  idx_t * options;
54  idx_t * objval;
56  idx_t * part;
57 };
58 
60 #define BALANCE_CC 1
61 #define BALANCE_CCM 2
63 #define BALANCE_CC_O(c) c+1
65 
72 template<typename Graph>
73 class Metis
74 {
77 
79  Graph & g;
80 
82  size_t n_dec;
83 
85  real_t dist_tol = 1.05;
86 
92  void constructAdjList(Graph & g)
93  {
94  // create xadj and adjlist
95 
96  // create xadj, adjlist, vwgt, adjwgt and vsize
97  Mg.xadj = new idx_t[g.getNVertex() + 1];
98  Mg.adjncy = new idx_t[g.getNEdge()];
99 
101  size_t prev = 0;
102 
103  // actual position
104  size_t id = 0;
105 
106  // for each vertex calculate the position of the starting point in the adjacency list
107  for (size_t i = 0; i < g.getNVertex(); i++)
108  {
109  // Calculate the starting point in the adjacency list
110  Mg.xadj[id] = prev;
111 
112  // Create the adjacency list
113  for (size_t s = 0; s < g.getNChilds(i); s++)
114  {
115  Mg.adjncy[prev + s] = g.getChild(i, s);
116  }
117 
118  // update the position for the next vertex
119  prev += g.getNChilds(i);
120 
121  id++;
122  }
123 
124  // Fill the last point
125  Mg.xadj[id] = prev;
126  }
127 
134  {
135  // create xadj, adjlist, vwgt, adjwgt and vsize
136  Mg.xadj = new idx_t[g.getNVertex() + 1];
137  Mg.adjncy = new idx_t[g.getNEdge()];
138  Mg.vwgt = new idx_t[g.getNVertex()];
139  Mg.adjwgt = new idx_t[g.getNEdge()];
140  Mg.vsize = new idx_t[g.getNVertex()];
141 
143  size_t prev = 0;
144 
145  // actual position
146  size_t id = 0;
147 
148  // for each vertex calculate the position of the starting point in the adjacency list
149  for (size_t i = 0; i < g.getNVertex(); i++)
150  {
151  // Add weight to vertex and migration cost
152  Mg.vwgt[i] = g.vertex(i).template get<nm_v::computation>();
153  Mg.vwgt[i] = (Mg.vwgt[i] == 0)?1:Mg.vwgt[i];
154  Mg.vsize[i] = g.vertex(i).template get<nm_v::migration>();
155  Mg.vsize[i] = (Mg.vsize[i] == 0)?1:Mg.vsize[i];
156 
157  // Calculate the starting point in the adjacency list
158  Mg.xadj[id] = prev;
159 
160  // Create the adjacency list
161  for (size_t s = 0; s < g.getNChilds(i); s++)
162  {
163  Mg.adjncy[prev + s] = g.getChild(i, s);
164 
165  // zero values on Metis are dangerous
166  Mg.adjwgt[prev + s] = g.getChildEdge(i, s).template get<nm_e::communication>();
167  Mg.adjwgt[prev + s] = (Mg.adjwgt[prev + s] == 0)?1:Mg.adjwgt[prev + s];
168  }
169 
170  // update the position for the next vertex
171  prev += g.getNChilds(i);
172 
173  id++;
174  }
175 
176  // Fill the last point
177  Mg.xadj[id] = prev;
178  }
179 
180 public:
181 
191  Metis(Graph & g, size_t nc, bool useWeights)
192  :g(g),n_dec(0)
193  {
194  initMetisGraph(nc,useWeights);
195  }
196 
205  Metis(Graph & g, size_t nc)
206  :g(g),n_dec(0)
207  {
208  initMetisGraph(nc,false);
209  }
210 
219  Metis(Graph & g)
220  :g(g),n_dec(0)
221  {
222  Mg.nvtxs = NULL;
223  Mg.ncon = NULL;
224  Mg.xadj = NULL;
225  Mg.adjncy = NULL;
226  Mg.vwgt = NULL;
227  Mg.adjwgt = NULL;
228  Mg.nparts = NULL;
229  Mg.tpwgts = NULL;
230  Mg.ubvec = NULL;
231  Mg.options = NULL;
232  Mg.objval = NULL;
233  Mg.part = NULL;
234  Mg.vsize = NULL;
235  }
236 
237 
244  void initMetisGraph(int nc, bool useWeights)
245  {
246 
247  // Get the number of vertex
248 
249  Mg.nvtxs = new idx_t[1];
250  Mg.nvtxs[0] = g.getNVertex();
251 
252  // Set the number of constrains
253 
254  Mg.ncon = new idx_t[1];
255  Mg.ncon[0] = 1;
256 
257  // Set to null the weight of the vertex
258 
259  Mg.vwgt = NULL;
260 
261  // Put the total communication size to NULL
262 
263  Mg.vsize = NULL;
264 
265  // Set to null the weight of the edge
266 
267  Mg.adjwgt = NULL;
268 
269  // construct the adjacency list
270  if (useWeights)
272  else
274 
275  // Set the total number of partitions
276 
277  Mg.nparts = new idx_t[1];
278  Mg.nparts[0] = nc;
279 
280  // Set to null the desired weight for each partition (one for each constrain)
281 
282  Mg.tpwgts = NULL;
283 
285 
286  Mg.ubvec = NULL;
287 
289 
290  Mg.options = NULL;
291 
293  Mg.objval = new idx_t[1];
294 
296  Mg.part = new idx_t[g.getNVertex()];
297 
298  for (size_t i = 0; i < g.getNVertex(); i++)
299  Mg.part[i] = 0;
300  }
301 
308  {
309  // Deallocate the Mg structure
310  if (Mg.nvtxs != NULL)
311  {
312  delete[] Mg.nvtxs;
313  }
314 
315  if (Mg.ncon != NULL)
316  {
317  delete[] Mg.ncon;
318  }
319 
320  if (Mg.xadj != NULL)
321  {
322  delete[] Mg.xadj;
323  }
324 
325  if (Mg.adjncy != NULL)
326  {
327  delete[] Mg.adjncy;
328  }
329  if (Mg.vwgt != NULL)
330  {
331  delete[] Mg.vwgt;
332  }
333  if (Mg.adjwgt != NULL)
334  {
335  delete[] Mg.adjwgt;
336  }
337  if (Mg.nparts != NULL)
338  {
339  delete[] Mg.nparts;
340  }
341  if (Mg.tpwgts != NULL)
342  {
343  delete[] Mg.tpwgts;
344  }
345  if (Mg.ubvec != NULL)
346  {
347  delete[] Mg.ubvec;
348  }
349  if (Mg.options != NULL)
350  {
351  delete[] Mg.options;
352  }
353  if (Mg.objval != NULL)
354  {
355  delete[] Mg.objval;
356  }
357  if (Mg.part != NULL)
358  {
359  delete[] Mg.part;
360  }
361  }
362 
369  template<unsigned int i>
370  void decompose()
371  {
372  if (Mg.nparts[0] != 1)
373  {
374  // Decompose
375  METIS_PartGraphRecursive(Mg.nvtxs, Mg.ncon, Mg.xadj, Mg.adjncy, Mg.vwgt, Mg.vsize, Mg.adjwgt, Mg.nparts, Mg.tpwgts, Mg.ubvec, Mg.options, Mg.objval, Mg.part);
376 
377  // vertex id
378 
379  size_t id = 0;
380 
381  // For each vertex store the processor that contain the data
382 
383  auto it = g.getVertexIterator();
384 
385  while (it.isNext())
386  {
387  g.vertex(it).template get<i>() = Mg.part[id];
388 
389  ++id;
390  ++it;
391  }
392  }
393  else
394  {
395  // Trivially assign all the domains to the processor 0
396 
397  auto it = g.getVertexIterator();
398 
399  while (it.isNext())
400  {
401  g.vertex(it).template get<i>() = 0;
402 
403  ++it;
404  }
405  }
406 
407  n_dec++;
408  }
409 
416  template<unsigned int i, typename Graph_part>
417  void decompose(Graph_part & gp)
418  {
419  // Decompose
420  METIS_PartGraphRecursive(Mg.nvtxs, Mg.ncon, Mg.xadj, Mg.adjncy, Mg.vwgt, Mg.vsize, Mg.adjwgt, Mg.nparts, Mg.tpwgts, Mg.ubvec, Mg.options, Mg.objval, Mg.part);
421 
422  // vertex id
423 
424  size_t id = 0;
425 
426  // For each vertex store the processor that contain the data
427 
428  auto it = gp.getVertexIterator();
429 
430  while (it.isNext())
431  {
432  gp.vertex(it).template get<i>() = Mg.part[id];
433 
434  ++id;
435  ++it;
436  }
437 
438  n_dec++;
439  }
440 
449  void onTest(bool testing)
450  {
451  if (testing == false)
452  return;
453 
454  if (Mg.options == NULL)
455  {
456  // allocate
457  Mg.options = new idx_t[METIS_NOPTIONS];
458 
459  // set default options
460  METIS_SetDefaultOptions(Mg.options);
461  }
462 
463  Mg.options[METIS_OPTION_SEED] = 0;
464  }
465 
471  void setDistTol(real_t tol)
472  {
473  dist_tol = tol;
474  }
475 
481  size_t get_ndec()
482  {
483  return n_dec;
484  }
485 
490  void inc_dec()
491  {
492  n_dec++;
493  }
494 };
495 
496 #endif
real_t dist_tol
Distribution tolerance.
Definition: metis_util.hpp:85
idx_t * nvtxs
The number of vertices in the graph.
Definition: metis_util.hpp:24
void onTest(bool testing)
It set Metis on test.
Definition: metis_util.hpp:449
void decompose()
Decompose the graph.
Definition: metis_util.hpp:370
void constructAdjList(Graph &g)
Construct Adjacency list.
Definition: metis_util.hpp:92
Helper class to define Metis graph.
Definition: metis_util.hpp:73
idx_t * vwgt
Array that store the weight for each vertex.
Definition: metis_util.hpp:39
idx_t * nparts
number of part to partition the graph
Definition: metis_util.hpp:46
Metis(Graph &g, size_t nc)
Constructor.
Definition: metis_util.hpp:205
idx_t * ncon
Definition: metis_util.hpp:30
void constructAdjListWithWeights(Graph &g)
Construct Adjacency list.
Definition: metis_util.hpp:133
Metis(Graph &g, size_t nc, bool useWeights)
Constructor.
Definition: metis_util.hpp:191
Metis_graph Mg
Graph in metis reppresentation.
Definition: metis_util.hpp:76
void setDistTol(real_t tol)
Distribution tolerance.
Definition: metis_util.hpp:471
idx_t * xadj
For each vertex it store the adjacency lost start for the vertex i.
Definition: metis_util.hpp:33
~Metis()
destructor
Definition: metis_util.hpp:307
idx_t * objval
return the total comunication cost for each partition
Definition: metis_util.hpp:54
idx_t * adjncy
For each vertex it store a list of all neighborhood vertex.
Definition: metis_util.hpp:36
idx_t * part
Is a output vector containing the partition for each vertex.
Definition: metis_util.hpp:56
idx_t * options
Additional option for the graph partitioning.
Definition: metis_util.hpp:52
idx_t * adjwgt
The weight of the edge.
Definition: metis_util.hpp:44
Metis graph structure.
Definition: metis_util.hpp:21
size_t get_ndec()
Get the decomposition counter.
Definition: metis_util.hpp:481
Graph & g
Original graph.
Definition: metis_util.hpp:79
Metis(Graph &g)
Constructor.
Definition: metis_util.hpp:219
real_t * tpwgts
Desired weight for each partition (one for each constrain)
Definition: metis_util.hpp:48
void inc_dec()
Increment the decomposition counter.
Definition: metis_util.hpp:490
void decompose(Graph_part &gp)
Decompose the graph.
Definition: metis_util.hpp:417
idx_t * vsize
Array of the vertex size, basically is the total communication amount.
Definition: metis_util.hpp:42
real_t * ubvec
For each partition load imbalance tollerated.
Definition: metis_util.hpp:50
void initMetisGraph(int nc, bool useWeights)
Initialize the METIS graph.
Definition: metis_util.hpp:244
size_t n_dec
indicate how many time decompose/refine/re-decompose has been called
Definition: metis_util.hpp:82