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