Zoltan2
Zoltan2_GraphMetricsUtility.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 
49 #ifndef ZOLTAN2_GRAPHICMETRICVALUESUTILITY_HPP
50 #define ZOLTAN2_GRAPHICMETRICVALUESUTILITY_HPP
51 
54 #include <zoltan_dd.h>
55 #include <Zoltan2_TPLTraits.hpp>
56 #include <map>
57 #include <vector>
58 
59 namespace Zoltan2{
60 
61 
62 template <typename Adapter, typename MachineRep>
64  const RCP<const Environment> &env,
65  const RCP<const Comm<int> > &comm,
66  const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graph,
67  const ArrayView<const typename Adapter::part_t> &parts,
68  typename Adapter::part_t &numParts,
69  ArrayRCP<RCP<BaseClassMetrics<typename Adapter::scalar_t> > > &metrics,
70  ArrayRCP<typename Adapter::scalar_t> &globalSums,
71  const RCP <const MachineRep> machine)
72 {
73  env->debug(DETAILED_STATUS, "Entering globalWeightedCutsMessagesHopsByPart");
75  // Initialize return values
76 
77  typedef typename Adapter::lno_t t_lno_t;
78  typedef typename Adapter::gno_t t_gno_t;
79  typedef typename Adapter::offset_t t_offset_t;
80  typedef typename Adapter::scalar_t t_scalar_t;
81  typedef typename Adapter::part_t part_t;
82  typedef typename Adapter::node_t t_node_t;
83 
84 
86 
87  t_lno_t localNumVertices = graph->getLocalNumVertices();
88  t_lno_t localNumEdges = graph->getLocalNumEdges();
89 
90  ArrayView<const t_gno_t> Ids;
91  ArrayView<t_input_t> v_wghts;
92  graph->getVertexList(Ids, v_wghts);
93 
94  typedef GraphMetrics<t_scalar_t> gm_t;
95 
96  //get the edge ids, and weights
97  ArrayView<const t_gno_t> edgeIds;
98  ArrayView<const t_offset_t> offsets;
99  ArrayView<t_input_t> e_wgts;
100  graph->getEdgeList(edgeIds, offsets, e_wgts);
101 
102 
103  std::vector <t_scalar_t> edge_weights;
104  int numWeightPerEdge = graph->getNumWeightsPerEdge();
105 
106  int numMetrics = 4; // "edge cuts", messages, hops, weighted hops
107  if (numWeightPerEdge) numMetrics += numWeightPerEdge * 2; // "weight n", weighted hops per weight n
108 
109  // add some more metrics to the array
110  auto next = metrics.size(); // where we begin filling
111  for (auto n = 0; n < numMetrics; ++n) {
112  addNewMetric<gm_t, t_scalar_t>(env, metrics);
113  }
114 
115  std::vector <part_t> e_parts (localNumEdges);
116 #ifdef HAVE_ZOLTAN2_MPI
117  if (comm->getSize() > 1)
118  {
119  Zoltan_DD_Struct *dd = NULL;
120 
121  MPI_Comm mpicomm = Teuchos::getRawMpiComm(*comm);
123 
124  int debug_level = 0;
125  Zoltan_DD_Create(&dd, mpicomm,
126  size_gnot, 0,
127  sizeof(part_t), localNumVertices, debug_level);
128 
129  ZOLTAN_ID_PTR ddnotneeded = NULL; // Local IDs not needed
130  Zoltan_DD_Update(
131  dd,
132  (localNumVertices ? (ZOLTAN_ID_PTR) Ids.getRawPtr() : NULL),
133  ddnotneeded,
134  (localNumVertices ? (char *) &(parts[0]) : NULL),
135  NULL,
136  int(localNumVertices));
137 
138  Zoltan_DD_Find(
139  dd,
140  (localNumEdges ? (ZOLTAN_ID_PTR) edgeIds.getRawPtr() : NULL),
141  ddnotneeded,
142  (localNumEdges ? (char *)&(e_parts[0]) : NULL),
143  NULL,
144  localNumEdges,
145  NULL
146  );
147  Zoltan_DD_Destroy(&dd);
148  } else
149 #endif
150  {
151 
152  std::map<t_gno_t,t_lno_t> global_id_to_local_index;
153 
154  //else everything is local.
155  //we need a globalid to local index conversion.
156  //this does not exists till this point, so we need to create one.
157  for (t_lno_t i = 0; i < localNumVertices; ++i){
158  //at the local index i, we have the global index Ids[i].
159  //so write i, to Ids[i] index of the vector.
160  global_id_to_local_index[Ids[i]] = i;
161  }
162 
163  for (t_lno_t i = 0; i < localNumEdges; ++i){
164  t_gno_t ei = edgeIds[i];
165  //ei is the global index of the neighbor one.
166  part_t p = parts[global_id_to_local_index[ei]];
167  e_parts[i] = p;
168  }
169  }
170 
171  RCP<const Teuchos::Comm<int> > tcomm = comm;
172 
173  env->timerStart(MACRO_TIMERS, "Communication Graph Create");
174  {
175  //get the vertices in each part in my part.
176  std::vector <t_lno_t> part_begins(numParts, -1);
177  std::vector <t_lno_t> part_nexts(localNumVertices, -1);
178 
179  //cluster vertices according to their parts.
180  //create local part graph.
181  for (t_lno_t i = 0; i < localNumVertices; ++i){
182  part_t ap = parts[i];
183  part_nexts[i] = part_begins[ap];
184  part_begins[ap] = i;
185  }
186 
187 
188  for (int weight_index = -1; weight_index < numWeightPerEdge ; ++weight_index){
189 
190  //MD: these two should be part_t.
191  //but we dont want to compile tpetra from the beginning.
192  //This can be changed when directory is updated.
193  typedef t_lno_t local_part_type;
194  typedef t_gno_t global_part_type;
195 
196  typedef Tpetra::Map<local_part_type, global_part_type, t_node_t> map_t;
197  Teuchos::RCP<const map_t> map = Teuchos::rcp (new map_t (numParts, 0, tcomm));
198 
199  typedef Tpetra::CrsMatrix<t_scalar_t, local_part_type, global_part_type, t_node_t> tcrsMatrix_t;
200  Teuchos::RCP<tcrsMatrix_t> tMatrix(new tcrsMatrix_t (map, 0));
201 
202 
203  std::vector <global_part_type> part_neighbors (numParts);
204 
205  std::vector <t_scalar_t> part_neighbor_weights(numParts, 0);
206  std::vector <t_scalar_t> part_neighbor_weights_ordered(numParts);
207 
208  //coarsen for all vertices in my part in order with parts.
209  for (global_part_type i = 0; i < (global_part_type) numParts; ++i){
210  part_t num_neighbor_parts = 0;
211  t_lno_t v = part_begins[i];
212  //get part i, and first vertex in this part v.
213  while (v != -1){
214  //now get the neightbors of v.
215  for (t_offset_t j = offsets[v]; j < offsets[v+1]; ++j){
216  //get the part of the second vertex.
217  part_t ep = e_parts[j];
218 
219  t_scalar_t ew = 1;
220  if (weight_index > -1){
221  ew = e_wgts[weight_index][j];
222  }
223  //add it to my local part neighbors for part i.
224  if (part_neighbor_weights[ep] < 0.00001){
225  part_neighbors[num_neighbor_parts++] = ep;
226  }
227  part_neighbor_weights[ep] += ew;
228  }
229  v = part_nexts[v];
230  }
231 
232  //now get the part list.
233  for (t_lno_t j = 0; j < num_neighbor_parts; ++j){
234  part_t neighbor_part = part_neighbors[j];
235  part_neighbor_weights_ordered[j] = part_neighbor_weights[neighbor_part];
236  part_neighbor_weights[neighbor_part] = 0;
237  }
238 
239  //insert it to tpetra crsmatrix.
240  if (num_neighbor_parts > 0){
241  Teuchos::ArrayView<const global_part_type> destinations(&(part_neighbors[0]), num_neighbor_parts);
242  Teuchos::ArrayView<const t_scalar_t> vals(&(part_neighbor_weights_ordered[0]), num_neighbor_parts);
243  tMatrix->insertGlobalValues (i,destinations, vals);
244  }
245  }
246  tMatrix->fillComplete ();
247  local_part_type num_local_parts = map->getNodeNumElements();
248 
249  Array<global_part_type> Indices;
250  Array<t_scalar_t> Values;
251 
252  t_scalar_t max_edge_cut = 0;
253  t_scalar_t total_edge_cut = 0;
254  global_part_type max_message = 0;
255  global_part_type total_message = 0;
256 
257  global_part_type total_hop_count = 0;
258  t_scalar_t total_weighted_hop_count = 0;
259  global_part_type max_hop_count = 0;
260  t_scalar_t max_weighted_hop_count = 0;
261 
262  for (local_part_type i=0; i < num_local_parts; i++) {
263 
264  const global_part_type globalRow = map->getGlobalElement(i);
265  size_t NumEntries = tMatrix->getNumEntriesInGlobalRow (globalRow);
266  Indices.resize (NumEntries);
267  Values.resize (NumEntries);
268  tMatrix->getGlobalRowCopy (globalRow,Indices(),Values(),NumEntries);
269 
270  t_scalar_t part_edge_cut = 0;
271  global_part_type part_messages = 0;
272 
273  for (size_t j=0; j < NumEntries; j++){
274  if (Indices[j] != globalRow){
275  part_edge_cut += Values[j];
276  part_messages += 1;
277 
278  typename MachineRep::machine_pcoord_t hop_count = 0;
279  machine->getHopCount(globalRow, Indices[j], hop_count);
280 
281  global_part_type hop_counts = hop_count;
282  t_scalar_t weighted_hop_counts = hop_count * Values[j];
283 
284  total_hop_count += hop_counts;
285  total_weighted_hop_count += weighted_hop_counts;
286 
287  if (hop_counts > max_hop_count ){
288  max_hop_count = hop_counts;
289  }
290  if (weighted_hop_counts > max_weighted_hop_count ){
291  max_weighted_hop_count = weighted_hop_counts;
292  }
293  }
294  }
295  if (part_edge_cut > max_edge_cut){
296  max_edge_cut = part_edge_cut;
297  }
298  total_edge_cut += part_edge_cut;
299 
300  if (part_messages > max_message){
301  max_message = part_messages;
302  }
303  total_message += part_messages;
304 
305  }
306  t_scalar_t g_max_edge_cut = 0;
307  t_scalar_t g_total_edge_cut = 0;
308  global_part_type g_max_message = 0;
309  global_part_type g_total_message = 0;
310 
311 
312 
313  global_part_type g_total_hop_count = 0;
314  t_scalar_t g_total_weighted_hop_count = 0;
315  global_part_type g_max_hop_count = 0;
316  t_scalar_t g_max_weighted_hop_count = 0;
317 
318  try{
319 
320  Teuchos::reduceAll<int, t_scalar_t>(*comm,Teuchos::REDUCE_MAX,1,&max_edge_cut,&g_max_edge_cut);
321  Teuchos::reduceAll<int, global_part_type>(*comm,Teuchos::REDUCE_MAX,1,&max_message,&g_max_message);
322 
323  Teuchos::reduceAll<int, global_part_type>(*comm,Teuchos::REDUCE_MAX,1,&max_hop_count,&g_max_hop_count);
324  Teuchos::reduceAll<int, t_scalar_t>(*comm,Teuchos::REDUCE_MAX,1,&max_weighted_hop_count,&g_max_weighted_hop_count);
325 
326  Teuchos::reduceAll<int, t_scalar_t>(*comm,Teuchos::REDUCE_SUM,1,&total_edge_cut,&g_total_edge_cut);
327  Teuchos::reduceAll<int, global_part_type>(*comm,Teuchos::REDUCE_SUM,1,&total_message,&g_total_message);
328 
329  Teuchos::reduceAll<int, global_part_type>(*comm,Teuchos::REDUCE_SUM,1,&total_hop_count,&g_total_hop_count);
330  Teuchos::reduceAll<int, t_scalar_t>(*comm,Teuchos::REDUCE_SUM,1,&total_weighted_hop_count,&g_total_weighted_hop_count);
331 
332  }
334 
335 
336  if (weight_index == -1){
337  metrics[next]->setName("md edge cuts");
338  }
339  else {
340  std::ostringstream oss;
341  oss << "md weighted edge cuts" << weight_index;
342  metrics[next]->setName( oss.str());
343  }
344 
345  metrics[next]->setMetricValue("global maximum", g_max_edge_cut);
346  metrics[next]->setMetricValue("global sum", g_total_edge_cut);
347  next++;
348 
349  if (weight_index == -1){
350  metrics[next]->setName("message");
351  metrics[next]->setMetricValue("global maximum", g_max_message);
352  metrics[next]->setMetricValue("global sum", g_total_message);
353  next++;
354  }
355 
356 
357  if (weight_index == -1){
358  metrics[next]->setName("hops (No Weight)");
359  metrics[next]->setMetricValue("global maximum", g_max_hop_count);
360  metrics[next]->setMetricValue("global sum", g_total_hop_count);
361  next++;
362  }
363 
364  std::ostringstream oss;
365  oss << "weighted hops" << weight_index;
366  metrics[next]->setName( oss.str());
367  metrics[next]->setMetricValue("global maximum", g_max_weighted_hop_count);
368  metrics[next]->setMetricValue("global sum", g_total_weighted_hop_count);
369  next++;
370 
371  }
372  }
373  env->timerStop(MACRO_TIMERS, "Communication Graph Create");
374 
375  env->debug(DETAILED_STATUS, "Exiting globalWeightedCutsMessagesHopsByPart");
376 }
377 
378 
379 template <typename Adapter>
381  const RCP<const Environment> &env,
382  const RCP<const Comm<int> > &comm,
383  const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graph,
384  const ArrayView<const typename Adapter::part_t> &parts,
385  typename Adapter::part_t &numParts,
386  ArrayRCP<RCP<BaseClassMetrics<typename Adapter::scalar_t> > > &metrics,
387  ArrayRCP<typename Adapter::scalar_t> &globalSums)
388 {
389  env->debug(DETAILED_STATUS, "Entering globalWeightedCutsMessagesByPart");
391  // Initialize return values
392 
393  typedef typename Adapter::lno_t t_lno_t;
394  typedef typename Adapter::gno_t t_gno_t;
395  typedef typename Adapter::offset_t t_offset_t;
396  typedef typename Adapter::scalar_t t_scalar_t;
397  typedef typename Adapter::part_t part_t;
398  typedef typename Adapter::node_t t_node_t;
399 
400 
402 
403  t_lno_t localNumVertices = graph->getLocalNumVertices();
404  t_lno_t localNumEdges = graph->getLocalNumEdges();
405 
406  ArrayView<const t_gno_t> Ids;
407  ArrayView<t_input_t> v_wghts;
408  graph->getVertexList(Ids, v_wghts);
409 
410  typedef GraphMetrics<t_scalar_t> gm_t;
411 
412  //get the edge ids, and weights
413  ArrayView<const t_gno_t> edgeIds;
414  ArrayView<const t_offset_t> offsets;
415  ArrayView<t_input_t> e_wgts;
416  graph->getEdgeList(edgeIds, offsets, e_wgts);
417 
418 
419  std::vector <t_scalar_t> edge_weights;
420  int numWeightPerEdge = graph->getNumWeightsPerEdge();
421 
422  int numMetrics = 2; // "edge cuts", messages
423  if (numWeightPerEdge) numMetrics += numWeightPerEdge; // "weight n"
424 
425  // add some more metrics to the array
426  auto next = metrics.size(); // where we begin filling
427  for (auto n = 0; n < numMetrics; ++n) {
428  addNewMetric<gm_t, t_scalar_t>(env, metrics);
429  }
430 
431  std::vector <part_t> e_parts (localNumEdges);
432 #ifdef HAVE_ZOLTAN2_MPI
433  if (comm->getSize() > 1)
434  {
435  Zoltan_DD_Struct *dd = NULL;
436 
437  MPI_Comm mpicomm = Teuchos::getRawMpiComm(*comm);
439 
440  int debug_level = 0;
441  Zoltan_DD_Create(&dd, mpicomm,
442  size_gnot, 0,
443  sizeof(part_t), localNumVertices, debug_level);
444 
445  ZOLTAN_ID_PTR ddnotneeded = NULL; // Local IDs not needed
446  Zoltan_DD_Update(
447  dd,
448  (localNumVertices ? (ZOLTAN_ID_PTR) Ids.getRawPtr() : NULL),
449  ddnotneeded,
450  (localNumVertices ? (char *) &(parts[0]) : NULL),
451  NULL,
452  int(localNumVertices));
453 
454  Zoltan_DD_Find(
455  dd,
456  (localNumEdges ? (ZOLTAN_ID_PTR) edgeIds.getRawPtr() : NULL),
457  ddnotneeded,
458  (localNumEdges ? (char *)&(e_parts[0]) : NULL),
459  NULL,
460  localNumEdges,
461  NULL
462  );
463  Zoltan_DD_Destroy(&dd);
464  } else
465 #endif
466  {
467 
468  std::map<t_gno_t,t_lno_t> global_id_to_local_index;
469 
470  //else everything is local.
471  //we need a globalid to local index conversion.
472  //this does not exists till this point, so we need to create one.
473  for (t_lno_t i = 0; i < localNumVertices; ++i){
474  //at the local index i, we have the global index Ids[i].
475  //so write i, to Ids[i] index of the vector.
476  global_id_to_local_index[Ids[i]] = i;
477  }
478 
479  for (t_lno_t i = 0; i < localNumEdges; ++i){
480  t_gno_t ei = edgeIds[i];
481  //ei is the global index of the neighbor one.
482  part_t p = parts[global_id_to_local_index[ei]];
483  e_parts[i] = p;
484  }
485  }
486 
487  RCP<const Teuchos::Comm<int> > tcomm = comm;
488 
489  env->timerStart(MACRO_TIMERS, "Communication Graph Create");
490  {
491  //get the vertices in each part in my part.
492  std::vector <t_lno_t> part_begins(numParts, -1);
493  std::vector <t_lno_t> part_nexts(localNumVertices, -1);
494 
495  //cluster vertices according to their parts.
496  //create local part graph.
497  for (t_lno_t i = 0; i < localNumVertices; ++i){
498  part_t ap = parts[i];
499  part_nexts[i] = part_begins[ap];
500  part_begins[ap] = i;
501  }
502 
503  for (int weight_index = -1; weight_index < numWeightPerEdge ; ++weight_index){
504 
505  //MD: these two should be part_t.
506  //but we dont want to compile tpetra from the beginning.
507  //This can be changed when directory is updated.
508  typedef t_lno_t local_part_type;
509  typedef t_gno_t global_part_type;
510 
511  typedef Tpetra::Map<local_part_type, global_part_type, t_node_t> map_t;
512  Teuchos::RCP<const map_t> map = Teuchos::rcp (new map_t (numParts, 0, tcomm));
513 
514  typedef Tpetra::CrsMatrix<t_scalar_t, local_part_type, global_part_type, t_node_t> tcrsMatrix_t;
515  Teuchos::RCP<tcrsMatrix_t> tMatrix(new tcrsMatrix_t (map, 0));
516 
517 
518  std::vector <global_part_type> part_neighbors (numParts);
519 
520  std::vector <t_scalar_t> part_neighbor_weights(numParts, 0);
521  std::vector <t_scalar_t> part_neighbor_weights_ordered(numParts);
522 
523  //coarsen for all vertices in my part in order with parts.
524  for (global_part_type i = 0; i < (global_part_type) numParts; ++i){
525  part_t num_neighbor_parts = 0;
526  t_lno_t v = part_begins[i];
527  //get part i, and first vertex in this part v.
528  while (v != -1){
529  //now get the neightbors of v.
530  for (t_offset_t j = offsets[v]; j < offsets[v+1]; ++j){
531  //get the part of the second vertex.
532  part_t ep = e_parts[j];
533 
534  t_scalar_t ew = 1;
535  if (weight_index > -1){
536  ew = e_wgts[weight_index][j];
537  }
538  //add it to my local part neighbors for part i.
539  if (part_neighbor_weights[ep] < 0.00001){
540  part_neighbors[num_neighbor_parts++] = ep;
541  }
542  part_neighbor_weights[ep] += ew;
543  }
544  v = part_nexts[v];
545  }
546 
547  //now get the part list.
548  for (t_lno_t j = 0; j < num_neighbor_parts; ++j){
549  part_t neighbor_part = part_neighbors[j];
550  part_neighbor_weights_ordered[j] = part_neighbor_weights[neighbor_part];
551  part_neighbor_weights[neighbor_part] = 0;
552  }
553 
554  //insert it to tpetra crsmatrix.
555  if (num_neighbor_parts > 0){
556  Teuchos::ArrayView<const global_part_type> destinations(&(part_neighbors[0]), num_neighbor_parts);
557  Teuchos::ArrayView<const t_scalar_t> vals(&(part_neighbor_weights_ordered[0]), num_neighbor_parts);
558  tMatrix->insertGlobalValues (i,destinations, vals);
559  }
560  }
561  tMatrix->fillComplete ();
562  local_part_type num_local_parts = map->getNodeNumElements();
563 
564  Array<global_part_type> Indices;
565  Array<t_scalar_t> Values;
566 
567  t_scalar_t max_edge_cut = 0;
568  t_scalar_t total_edge_cut = 0;
569  global_part_type max_message = 0;
570  global_part_type total_message = 0;
571 
572  for (local_part_type i=0; i < num_local_parts; i++) {
573 
574  const global_part_type globalRow = map->getGlobalElement(i);
575  size_t NumEntries = tMatrix->getNumEntriesInGlobalRow (globalRow);
576  Indices.resize (NumEntries);
577  Values.resize (NumEntries);
578  tMatrix->getGlobalRowCopy (globalRow,Indices(),Values(),NumEntries);
579 
580  t_scalar_t part_edge_cut = 0;
581  global_part_type part_messages = 0;
582  for (size_t j=0; j < NumEntries; j++){
583  if (Indices[j] != globalRow){
584  part_edge_cut += Values[j];
585  part_messages += 1;
586  }
587  }
588  if (part_edge_cut > max_edge_cut){
589  max_edge_cut = part_edge_cut;
590  }
591  total_edge_cut += part_edge_cut;
592 
593  if (part_messages > max_message){
594  max_message = part_messages;
595  }
596  total_message += part_messages;
597 
598  }
599  t_scalar_t g_max_edge_cut = 0;
600  t_scalar_t g_total_edge_cut = 0;
601  global_part_type g_max_message = 0;
602  global_part_type g_total_message = 0;
603 
604  try{
605 
606  Teuchos::reduceAll<int, t_scalar_t>(*comm,Teuchos::REDUCE_MAX,1,&max_edge_cut,&g_max_edge_cut);
607  Teuchos::reduceAll<int, global_part_type>(*comm,Teuchos::REDUCE_MAX,1,&max_message,&g_max_message);
608 
609  Teuchos::reduceAll<int, t_scalar_t>(*comm,Teuchos::REDUCE_SUM,1,&total_edge_cut,&g_total_edge_cut);
610  Teuchos::reduceAll<int, global_part_type>(*comm,Teuchos::REDUCE_SUM,1,&total_message,&g_total_message);
611  }
613 
614  if (weight_index == -1){
615  metrics[next]->setName("md edge cuts");
616  }
617  else {
618  std::ostringstream oss;
619  oss << "md weight " << weight_index;
620  metrics[next]->setName( oss.str());
621  }
622  metrics[next]->setMetricValue("global maximum", g_max_edge_cut);
623  metrics[next]->setMetricValue("global sum", g_total_edge_cut);
624  next++;
625  if (weight_index == -1){
626  metrics[next]->setName("md message");
627  metrics[next]->setMetricValue("global maximum", g_max_message);
628  metrics[next]->setMetricValue("global sum", g_total_message);
629  next++;
630  }
631  }
632  }
633  env->timerStop(MACRO_TIMERS, "Communication Graph Create");
634 
635  env->debug(DETAILED_STATUS, "Exiting globalWeightedCutsMessagesByPart");
636 }
637 
666 template <typename Adapter>
668  const RCP<const Environment> &env,
669  const RCP<const Comm<int> > &comm,
670  const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graph,
671  const ArrayView<const typename Adapter::part_t> &part,
672  typename Adapter::part_t &numParts,
673  ArrayRCP<RCP<BaseClassMetrics<typename Adapter::scalar_t> > > &metrics,
674  ArrayRCP<typename Adapter::scalar_t> &globalSums)
675 {
676  env->debug(DETAILED_STATUS, "Entering globalWeightedCutsByPart");
678  // Initialize return values
679 
680  numParts = 0;
681 
682  int ewgtDim = graph->getNumWeightsPerEdge();
683 
684  int numMetrics = 1; // "edge cuts"
685  if (ewgtDim) numMetrics += ewgtDim; // "weight n"
686 
687  typedef typename Adapter::scalar_t scalar_t;
688  typedef typename Adapter::gno_t gno_t;
689  typedef typename Adapter::lno_t lno_t;
690  typedef typename Adapter::offset_t offset_t;
691  typedef typename Adapter::node_t node_t;
692  typedef typename Adapter::part_t part_t;
693  typedef StridedData<lno_t, scalar_t> input_t;
694 
695  typedef Tpetra::CrsMatrix<part_t,lno_t,gno_t,node_t> sparse_matrix_type;
696  typedef Tpetra::Vector<part_t,lno_t,gno_t,node_t> vector_t;
697  typedef Tpetra::Map<lno_t, gno_t, node_t> map_type;
698  typedef Tpetra::global_size_t GST;
699  const GST INVALID = Teuchos::OrdinalTraits<GST>::invalid ();
700 
701  using Teuchos::as;
702 
703  auto next = metrics.size(); // where we begin filling
704  typedef GraphMetrics<scalar_t> gm_t;
705  for (auto n = 0; n < numMetrics; ++n) {
706  addNewMetric<gm_t, scalar_t>(env, metrics);
707  }
708 
710  // Figure out the global number of parts in use.
711  // Verify number of vertex weights is the same everywhere.
712 
713  lno_t localNumObj = part.size();
714  part_t localNum[2], globalNum[2];
715  localNum[0] = static_cast<part_t>(ewgtDim);
716  localNum[1] = 0;
717 
718  for (lno_t i=0; i < localNumObj; i++)
719  if (part[i] > localNum[1]) localNum[1] = part[i];
720 
721  try{
722  reduceAll<int, part_t>(*comm, Teuchos::REDUCE_MAX, 2,
723  localNum, globalNum);
724  }
726 
727  env->globalBugAssertion(__FILE__,__LINE__,
728  "inconsistent number of edge weights",
729  globalNum[0] == localNum[0], DEBUG_MODE_ASSERTION, comm);
730 
731  part_t nparts = globalNum[1] + 1;
732 
733  part_t globalSumSize = nparts * numMetrics;
734  scalar_t * sumBuf = new scalar_t [globalSumSize];
735  env->localMemoryAssertion(__FILE__, __LINE__, globalSumSize, sumBuf);
736  globalSums = arcp(sumBuf, 0, globalSumSize);
737 
739  // Calculate the local totals by part.
740 
741  scalar_t *localBuf = new scalar_t [globalSumSize];
742  env->localMemoryAssertion(__FILE__,__LINE__,globalSumSize,localBuf);
743  memset(localBuf, 0, sizeof(scalar_t) * globalSumSize);
744 
745  scalar_t *cut = localBuf; // # of cuts
746 
747  ArrayView<const gno_t> Ids;
748  ArrayView<input_t> vwgts;
749  //size_t nv =
750  graph->getVertexList(Ids, vwgts);
751 
752  ArrayView<const gno_t> edgeIds;
753  ArrayView<const offset_t> offsets;
754  ArrayView<input_t> wgts;
755  //size_t numLocalEdges =
756  graph->getEdgeList(edgeIds, offsets, wgts);
757  // **************************************************************************
758  // *************************** BUILD MAP FOR ADJS ***************************
759  // **************************************************************************
760 
761  RCP<const map_type> vertexMapG;
762 
763  // Build a list of the global vertex ids...
764  gno_t min = std::numeric_limits<gno_t>::max();
765  offset_t maxcols = 0;
766  for (lno_t i = 0; i < localNumObj; ++i) {
767  if (Ids[i] < min) min = Ids[i];
768  offset_t ncols = offsets[i+1] - offsets[i];
769  if (ncols > maxcols) maxcols = ncols;
770  }
771 
772  gno_t gmin;
773  Teuchos::reduceAll<int, gno_t>(*comm,Teuchos::REDUCE_MIN,1,&min,&gmin);
774 
775  //Generate Map for vertex
776  vertexMapG = rcp(new map_type(INVALID, Ids, gmin, comm));
777 
778  // **************************************************************************
779  // ************************** BUILD GRAPH FOR ADJS **************************
780  // **************************************************************************
781 
782  //MD:Zoltan Directory could be used instead of adjMatrix.
783 
784  RCP<sparse_matrix_type> adjsMatrix;
785 
786  // Construct Tpetra::CrsGraph objects.
787  adjsMatrix = rcp (new sparse_matrix_type (vertexMapG, 0));
788 
789  Array<part_t> justOneA(maxcols, 1);
790 
791  for (lno_t localElement=0; localElement<localNumObj; ++localElement){
792  // Insert all columns for global row Ids[localElement]
793  offset_t ncols = offsets[localElement+1] - offsets[localElement];
794  adjsMatrix->insertGlobalValues(Ids[localElement],
795  edgeIds(offsets[localElement], ncols),
796  justOneA(0, ncols));
797  }
798 
799  //Fill-complete adjs Graph
800  adjsMatrix->fillComplete ();
801 
802  // Compute part
803  RCP<vector_t> scaleVec = Teuchos::rcp( new vector_t(vertexMapG,false) );
804  for (lno_t localElement=0; localElement<localNumObj; ++localElement) {
805  scaleVec->replaceLocalValue(localElement,part[localElement]);
806  }
807 
808  // Postmultiply adjsMatrix by part
809  adjsMatrix->rightScale(*scaleVec);
810  Array<gno_t> Indices;
811  Array<part_t> Values;
812 
813  for (lno_t i=0; i < localNumObj; i++) {
814  const gno_t globalRow = Ids[i];
815  size_t NumEntries = adjsMatrix->getNumEntriesInGlobalRow (globalRow);
816  Indices.resize (NumEntries);
817  Values.resize (NumEntries);
818  adjsMatrix->getGlobalRowCopy (globalRow,Indices(),Values(),NumEntries);
819 
820  for (size_t j=0; j < NumEntries; j++)
821  if (part[i] != Values[j])
822  cut[part[i]]++;
823  }
824 
825  if (numMetrics > 1) {
826 
827  scalar_t *wgt = localBuf + nparts; // weight 0
828 
829  // This code assumes the solution has the part ordered the
830  // same way as the user input. (Bug 5891 is resolved.)
831  for (int edim = 0; edim < ewgtDim; edim++){
832  for (lno_t i=0; i < localNumObj; i++) {
833  const gno_t globalRow = Ids[i];
834  size_t NumEntries = adjsMatrix->getNumEntriesInGlobalRow (globalRow);
835  Indices.resize (NumEntries);
836  Values.resize (NumEntries);
837  adjsMatrix->getGlobalRowCopy (globalRow,Indices(),Values(),NumEntries);
838 
839  for (size_t j=0; j < NumEntries; j++)
840  if (part[i] != Values[j])
841  wgt[part[i]] += wgts[edim][offsets[i] + j];
842  }
843  wgt += nparts; // individual weights
844  }
845  }
846 
848  // Obtain global totals by part.
849 
850  try{
851  reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_SUM, globalSumSize,
852  localBuf, sumBuf);
853  }
855 
856  delete [] localBuf;
857 
859  // Global max and sum over all parts
860 
861  cut = sumBuf; // # of cuts
862  scalar_t max=0, sum=0;
863 
864  ArrayView<scalar_t> cutVec(cut, nparts);
865  getStridedStats<scalar_t>(cutVec, 1, 0, max, sum);
866 
867  metrics[next]->setName("edge cuts");
868  metrics[next]->setMetricValue("global maximum", max);
869  metrics[next]->setMetricValue("global sum", sum);
870 
871  next++;
872 
873  if (numMetrics > 1){
874  scalar_t *wgt = sumBuf + nparts; // weight 0
875 
876  for (int edim=0; edim < ewgtDim; edim++){
877  ArrayView<scalar_t> fromVec(wgt, nparts);
878  getStridedStats<scalar_t>(fromVec, 1, 0, max, sum);
879 
880  std::ostringstream oss;
881  oss << "weight " << edim;
882 
883  metrics[next]->setName(oss.str());
884  metrics[next]->setMetricValue("global maximum", max);
885  metrics[next]->setMetricValue("global sum", sum);
886 
887  next++;
888  wgt += nparts; // individual weights
889  }
890  }
891 
892  numParts = nparts;
893 
894  env->debug(DETAILED_STATUS, "Exiting globalWeightedCutsByPart");
895 }
896 
899 template <typename scalar_t, typename part_t>
900 void printGraphMetricsHeader(std::ostream &os, part_t targetNumParts, part_t numParts )
901 {
902  os << "Graph Metrics: (" << numParts << " parts)";
903  os << std::endl;
904  if (targetNumParts != numParts) {
905  os << "Target number of parts is: " << targetNumParts << std::endl;
906  }
908 }
909 
912 template <typename scalar_t, typename part_t>
913 void printGraphMetrics(std::ostream &os, part_t targetNumParts, part_t numParts, const ArrayView<RCP<BaseClassMetrics<scalar_t> > > &infoList)
914 {
915  printGraphMetricsHeader<scalar_t, part_t>(os, targetNumParts, numParts);
916  for (int i=0; i < infoList.size(); i++) {
917  if (infoList[i]->getName() != METRICS_UNSET_STRING) {
918  infoList[i]->printLine(os);
919  }
920  }
921  os << std::endl;
922 }
923 
926 template <typename scalar_t, typename part_t>
927 void printGraphMetrics(std::ostream &os, part_t targetNumParts, part_t numParts, RCP<BaseClassMetrics<scalar_t>> metricValue)
928 {
929  printGraphMetricsHeader<scalar_t, part_t>(os, targetNumParts, numParts);
930  metricValue->printLine(os);
931 }
932 
933 } //namespace Zoltan2
934 
935 #endif
Zoltan2::GraphMetrics
Definition: Zoltan2_GraphMetrics.hpp:62
Zoltan2::GraphMetrics::printHeader
static void printHeader(std::ostream &os)
Print a standard header.
Definition: Zoltan2_GraphMetrics.hpp:114
Zoltan2_MetricUtility.hpp
Zoltan2::printGraphMetricsHeader
void printGraphMetricsHeader(std::ostream &os, part_t targetNumParts, part_t numParts)
Print out header info for graph metrics.
Definition: Zoltan2_GraphMetricsUtility.hpp:900
Zoltan2_ImbalanceMetrics.hpp
Zoltan2::MACRO_TIMERS
Time an algorithm (or other entity) as a whole.
Definition: Zoltan2_Parameters.hpp:120
INVALID
#define INVALID(STR)
Definition: GeometricGenerator.hpp:206
Zoltan2::DETAILED_STATUS
sub-steps, each method's entry and exit
Definition: Zoltan2_Parameters.hpp:101
xml2dox.vals
dictionary vals
Definition: xml2dox.py:186
part_t
SparseMatrixAdapter_t::part_t part_t
Definition: partitioningTree.cpp:74
Zoltan2::GraphModel::getLocalNumVertices
size_t getLocalNumVertices() const
Returns the number vertices on this process.
Definition: Zoltan2_GraphModel.hpp:141
Zoltan2::globalWeightedCutsByPart
void globalWeightedCutsByPart(const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, const RCP< const GraphModel< typename Adapter::base_adapter_t > > &graph, const ArrayView< const typename Adapter::part_t > &part, typename Adapter::part_t &numParts, ArrayRCP< RCP< BaseClassMetrics< typename Adapter::scalar_t > > > &metrics, ArrayRCP< typename Adapter::scalar_t > &globalSums)
Given the local partitioning, compute the global weighted cuts in each part.
Definition: Zoltan2_GraphMetricsUtility.hpp:667
Zoltan2::global_size_t
Tpetra::global_size_t global_size_t
Definition: Zoltan2_Standards.hpp:119
Zoltan2::DEBUG_MODE_ASSERTION
checks for logic errors
Definition: Zoltan2_Parameters.hpp:85
Z2_THROW_OUTSIDE_ERROR
#define Z2_THROW_OUTSIDE_ERROR(env)
Throw an error returned from outside the Zoltan2 library.
Definition: Zoltan2_Exceptions.hpp:64
Zoltan2::globalWeightedCutsMessagesHopsByPart
void globalWeightedCutsMessagesHopsByPart(const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, const RCP< const GraphModel< typename Adapter::base_adapter_t > > &graph, const ArrayView< const typename Adapter::part_t > &parts, typename Adapter::part_t &numParts, ArrayRCP< RCP< BaseClassMetrics< typename Adapter::scalar_t > > > &metrics, ArrayRCP< typename Adapter::scalar_t > &globalSums, const RCP< const MachineRep > machine)
Definition: Zoltan2_GraphMetricsUtility.hpp:63
METRICS_UNSET_STRING
#define METRICS_UNSET_STRING
Definition: Zoltan2_BaseClassMetrics.hpp:56
Zoltan2::printGraphMetrics
void printGraphMetrics(std::ostream &os, part_t targetNumParts, part_t numParts, const ArrayView< RCP< BaseClassMetrics< scalar_t > > > &infoList)
Print out list of graph metrics.
Definition: Zoltan2_GraphMetricsUtility.hpp:913
Zoltan2::TPL_Traits
Definition: Zoltan2_TPLTraits.hpp:70
Zoltan2::BaseClassMetrics
Definition: Zoltan2_BaseClassMetrics.hpp:62
Zoltan2
Definition: Zoltan2_AlgSerialGreedy.hpp:56
Zoltan2_TPLTraits.hpp
Traits class to handle conversions between gno_t/lno_t and TPL data types (e.g., ParMETIS's idx_t,...
tcrsMatrix_t
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tcrsMatrix_t
Definition: GraphModel.cpp:82
Zoltan2::globalWeightedCutsMessagesByPart
void globalWeightedCutsMessagesByPart(const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, const RCP< const GraphModel< typename Adapter::base_adapter_t > > &graph, const ArrayView< const typename Adapter::part_t > &parts, typename Adapter::part_t &numParts, ArrayRCP< RCP< BaseClassMetrics< typename Adapter::scalar_t > > > &metrics, ArrayRCP< typename Adapter::scalar_t > &globalSums)
Definition: Zoltan2_GraphMetricsUtility.hpp:380
Zoltan2::GraphModel
GraphModel defines the interface required for graph models.
Definition: Zoltan2_GraphModel.hpp:80
Zoltan2::StridedData
The StridedData class manages lists of weights or coordinates.
Definition: Zoltan2_StridedData.hpp:76