49 #ifndef ZOLTAN2_GRAPHICMETRICVALUESUTILITY_HPP
50 #define ZOLTAN2_GRAPHICMETRICVALUESUTILITY_HPP
54 #include <zoltan_dd.h>
62 template <
typename Adapter,
typename MachineRep>
64 const RCP<const Environment> &env,
65 const RCP<
const Comm<int> > &comm,
67 const ArrayView<const typename Adapter::part_t> &parts,
70 ArrayRCP<typename Adapter::scalar_t> &globalSums,
71 const RCP <const MachineRep> machine)
73 env->debug(
DETAILED_STATUS,
"Entering globalWeightedCutsMessagesHopsByPart");
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;
82 typedef typename Adapter::node_t t_node_t;
88 t_lno_t localNumEdges = graph->getLocalNumEdges();
90 ArrayView<const t_gno_t> Ids;
91 ArrayView<t_input_t> v_wghts;
92 graph->getVertexList(Ids, v_wghts);
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);
103 std::vector <t_scalar_t> edge_weights;
104 int numWeightPerEdge = graph->getNumWeightsPerEdge();
107 if (numWeightPerEdge) numMetrics += numWeightPerEdge * 2;
110 auto next = metrics.size();
111 for (
auto n = 0; n < numMetrics; ++n) {
112 addNewMetric<gm_t, t_scalar_t>(env, metrics);
115 std::vector <part_t> e_parts (localNumEdges);
116 #ifdef HAVE_ZOLTAN2_MPI
117 if (comm->getSize() > 1)
119 Zoltan_DD_Struct *dd = NULL;
121 MPI_Comm mpicomm = Teuchos::getRawMpiComm(*comm);
125 Zoltan_DD_Create(&dd, mpicomm,
127 sizeof(
part_t), localNumVertices, debug_level);
129 ZOLTAN_ID_PTR ddnotneeded = NULL;
132 (localNumVertices ? (ZOLTAN_ID_PTR) Ids.getRawPtr() : NULL),
134 (localNumVertices ? (
char *) &(parts[0]) : NULL),
136 int(localNumVertices));
140 (localNumEdges ? (ZOLTAN_ID_PTR) edgeIds.getRawPtr() : NULL),
142 (localNumEdges ? (
char *)&(e_parts[0]) : NULL),
147 Zoltan_DD_Destroy(&dd);
152 std::map<t_gno_t,t_lno_t> global_id_to_local_index;
157 for (t_lno_t i = 0; i < localNumVertices; ++i){
160 global_id_to_local_index[Ids[i]] = i;
163 for (t_lno_t i = 0; i < localNumEdges; ++i){
164 t_gno_t ei = edgeIds[i];
166 part_t p = parts[global_id_to_local_index[ei]];
171 RCP<const Teuchos::Comm<int> > tcomm = comm;
173 env->timerStart(
MACRO_TIMERS,
"Communication Graph Create");
176 std::vector <t_lno_t> part_begins(numParts, -1);
177 std::vector <t_lno_t> part_nexts(localNumVertices, -1);
181 for (t_lno_t i = 0; i < localNumVertices; ++i){
183 part_nexts[i] = part_begins[ap];
188 for (
int weight_index = -1; weight_index < numWeightPerEdge ; ++weight_index){
193 typedef t_lno_t local_part_type;
194 typedef t_gno_t global_part_type;
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));
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));
203 std::vector <global_part_type> part_neighbors (numParts);
205 std::vector <t_scalar_t> part_neighbor_weights(numParts, 0);
206 std::vector <t_scalar_t> part_neighbor_weights_ordered(numParts);
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];
215 for (t_offset_t j = offsets[v]; j < offsets[v+1]; ++j){
220 if (weight_index > -1){
221 ew = e_wgts[weight_index][j];
224 if (part_neighbor_weights[ep] < 0.00001){
225 part_neighbors[num_neighbor_parts++] = ep;
227 part_neighbor_weights[ep] += ew;
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;
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);
246 tMatrix->fillComplete ();
247 local_part_type num_local_parts = map->getNodeNumElements();
249 Array<global_part_type> Indices;
250 Array<t_scalar_t> Values;
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;
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;
262 for (local_part_type i=0; i < num_local_parts; i++) {
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);
270 t_scalar_t part_edge_cut = 0;
271 global_part_type part_messages = 0;
273 for (
size_t j=0; j < NumEntries; j++){
274 if (Indices[j] != globalRow){
275 part_edge_cut += Values[j];
278 typename MachineRep::machine_pcoord_t hop_count = 0;
279 machine->getHopCount(globalRow, Indices[j], hop_count);
281 global_part_type hop_counts = hop_count;
282 t_scalar_t weighted_hop_counts = hop_count * Values[j];
284 total_hop_count += hop_counts;
285 total_weighted_hop_count += weighted_hop_counts;
287 if (hop_counts > max_hop_count ){
288 max_hop_count = hop_counts;
290 if (weighted_hop_counts > max_weighted_hop_count ){
291 max_weighted_hop_count = weighted_hop_counts;
295 if (part_edge_cut > max_edge_cut){
296 max_edge_cut = part_edge_cut;
298 total_edge_cut += part_edge_cut;
300 if (part_messages > max_message){
301 max_message = part_messages;
303 total_message += part_messages;
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;
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;
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);
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);
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);
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);
336 if (weight_index == -1){
337 metrics[next]->setName(
"md edge cuts");
340 std::ostringstream oss;
341 oss <<
"md weighted edge cuts" << weight_index;
342 metrics[next]->setName( oss.str());
345 metrics[next]->setMetricValue(
"global maximum", g_max_edge_cut);
346 metrics[next]->setMetricValue(
"global sum", g_total_edge_cut);
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);
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);
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);
373 env->timerStop(
MACRO_TIMERS,
"Communication Graph Create");
375 env->debug(
DETAILED_STATUS,
"Exiting globalWeightedCutsMessagesHopsByPart");
379 template <
typename Adapter>
381 const RCP<const Environment> &env,
382 const RCP<
const Comm<int> > &comm,
384 const ArrayView<const typename Adapter::part_t> &parts,
387 ArrayRCP<typename Adapter::scalar_t> &globalSums)
389 env->debug(
DETAILED_STATUS,
"Entering globalWeightedCutsMessagesByPart");
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;
398 typedef typename Adapter::node_t t_node_t;
404 t_lno_t localNumEdges = graph->getLocalNumEdges();
406 ArrayView<const t_gno_t> Ids;
407 ArrayView<t_input_t> v_wghts;
408 graph->getVertexList(Ids, v_wghts);
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);
419 std::vector <t_scalar_t> edge_weights;
420 int numWeightPerEdge = graph->getNumWeightsPerEdge();
423 if (numWeightPerEdge) numMetrics += numWeightPerEdge;
426 auto next = metrics.size();
427 for (
auto n = 0; n < numMetrics; ++n) {
428 addNewMetric<gm_t, t_scalar_t>(env, metrics);
431 std::vector <part_t> e_parts (localNumEdges);
432 #ifdef HAVE_ZOLTAN2_MPI
433 if (comm->getSize() > 1)
435 Zoltan_DD_Struct *dd = NULL;
437 MPI_Comm mpicomm = Teuchos::getRawMpiComm(*comm);
441 Zoltan_DD_Create(&dd, mpicomm,
443 sizeof(
part_t), localNumVertices, debug_level);
445 ZOLTAN_ID_PTR ddnotneeded = NULL;
448 (localNumVertices ? (ZOLTAN_ID_PTR) Ids.getRawPtr() : NULL),
450 (localNumVertices ? (
char *) &(parts[0]) : NULL),
452 int(localNumVertices));
456 (localNumEdges ? (ZOLTAN_ID_PTR) edgeIds.getRawPtr() : NULL),
458 (localNumEdges ? (
char *)&(e_parts[0]) : NULL),
463 Zoltan_DD_Destroy(&dd);
468 std::map<t_gno_t,t_lno_t> global_id_to_local_index;
473 for (t_lno_t i = 0; i < localNumVertices; ++i){
476 global_id_to_local_index[Ids[i]] = i;
479 for (t_lno_t i = 0; i < localNumEdges; ++i){
480 t_gno_t ei = edgeIds[i];
482 part_t p = parts[global_id_to_local_index[ei]];
487 RCP<const Teuchos::Comm<int> > tcomm = comm;
489 env->timerStart(
MACRO_TIMERS,
"Communication Graph Create");
492 std::vector <t_lno_t> part_begins(numParts, -1);
493 std::vector <t_lno_t> part_nexts(localNumVertices, -1);
497 for (t_lno_t i = 0; i < localNumVertices; ++i){
499 part_nexts[i] = part_begins[ap];
503 for (
int weight_index = -1; weight_index < numWeightPerEdge ; ++weight_index){
508 typedef t_lno_t local_part_type;
509 typedef t_gno_t global_part_type;
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));
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));
518 std::vector <global_part_type> part_neighbors (numParts);
520 std::vector <t_scalar_t> part_neighbor_weights(numParts, 0);
521 std::vector <t_scalar_t> part_neighbor_weights_ordered(numParts);
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];
530 for (t_offset_t j = offsets[v]; j < offsets[v+1]; ++j){
535 if (weight_index > -1){
536 ew = e_wgts[weight_index][j];
539 if (part_neighbor_weights[ep] < 0.00001){
540 part_neighbors[num_neighbor_parts++] = ep;
542 part_neighbor_weights[ep] += ew;
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;
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);
561 tMatrix->fillComplete ();
562 local_part_type num_local_parts = map->getNodeNumElements();
564 Array<global_part_type> Indices;
565 Array<t_scalar_t> Values;
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;
572 for (local_part_type i=0; i < num_local_parts; i++) {
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);
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];
588 if (part_edge_cut > max_edge_cut){
589 max_edge_cut = part_edge_cut;
591 total_edge_cut += part_edge_cut;
593 if (part_messages > max_message){
594 max_message = part_messages;
596 total_message += part_messages;
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;
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);
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);
614 if (weight_index == -1){
615 metrics[next]->setName(
"md edge cuts");
618 std::ostringstream oss;
619 oss <<
"md weight " << weight_index;
620 metrics[next]->setName( oss.str());
622 metrics[next]->setMetricValue(
"global maximum", g_max_edge_cut);
623 metrics[next]->setMetricValue(
"global sum", g_total_edge_cut);
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);
633 env->timerStop(
MACRO_TIMERS,
"Communication Graph Create");
635 env->debug(
DETAILED_STATUS,
"Exiting globalWeightedCutsMessagesByPart");
666 template <
typename Adapter>
668 const RCP<const Environment> &env,
669 const RCP<
const Comm<int> > &comm,
671 const ArrayView<const typename Adapter::part_t> &part,
674 ArrayRCP<typename Adapter::scalar_t> &globalSums)
682 int ewgtDim = graph->getNumWeightsPerEdge();
685 if (ewgtDim) numMetrics += ewgtDim;
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;
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;
699 const GST
INVALID = Teuchos::OrdinalTraits<GST>::invalid ();
703 auto next = metrics.size();
705 for (
auto n = 0; n < numMetrics; ++n) {
706 addNewMetric<gm_t, scalar_t>(env, metrics);
713 lno_t localNumObj = part.size();
714 part_t localNum[2], globalNum[2];
715 localNum[0] = static_cast<part_t>(ewgtDim);
718 for (lno_t i=0; i < localNumObj; i++)
719 if (part[i] > localNum[1]) localNum[1] = part[i];
722 reduceAll<int, part_t>(*comm, Teuchos::REDUCE_MAX, 2,
723 localNum, globalNum);
727 env->globalBugAssertion(__FILE__,__LINE__,
728 "inconsistent number of edge weights",
731 part_t nparts = globalNum[1] + 1;
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);
741 scalar_t *localBuf =
new scalar_t [globalSumSize];
742 env->localMemoryAssertion(__FILE__,__LINE__,globalSumSize,localBuf);
743 memset(localBuf, 0,
sizeof(scalar_t) * globalSumSize);
745 scalar_t *cut = localBuf;
747 ArrayView<const gno_t> Ids;
748 ArrayView<input_t> vwgts;
750 graph->getVertexList(Ids, vwgts);
752 ArrayView<const gno_t> edgeIds;
753 ArrayView<const offset_t> offsets;
754 ArrayView<input_t> wgts;
756 graph->getEdgeList(edgeIds, offsets, wgts);
761 RCP<const map_type> vertexMapG;
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;
773 Teuchos::reduceAll<int, gno_t>(*comm,Teuchos::REDUCE_MIN,1,&min,&gmin);
776 vertexMapG = rcp(
new map_type(
INVALID, Ids, gmin, comm));
784 RCP<sparse_matrix_type> adjsMatrix;
787 adjsMatrix = rcp (
new sparse_matrix_type (vertexMapG, 0));
789 Array<part_t> justOneA(maxcols, 1);
791 for (lno_t localElement=0; localElement<localNumObj; ++localElement){
793 offset_t ncols = offsets[localElement+1] - offsets[localElement];
794 adjsMatrix->insertGlobalValues(Ids[localElement],
795 edgeIds(offsets[localElement], ncols),
800 adjsMatrix->fillComplete ();
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]);
809 adjsMatrix->rightScale(*scaleVec);
810 Array<gno_t> Indices;
811 Array<part_t> Values;
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);
820 for (
size_t j=0; j < NumEntries; j++)
821 if (part[i] != Values[j])
825 if (numMetrics > 1) {
827 scalar_t *wgt = localBuf + nparts;
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);
839 for (
size_t j=0; j < NumEntries; j++)
840 if (part[i] != Values[j])
841 wgt[part[i]] += wgts[edim][offsets[i] + j];
851 reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_SUM, globalSumSize,
862 scalar_t max=0, sum=0;
864 ArrayView<scalar_t> cutVec(cut, nparts);
865 getStridedStats<scalar_t>(cutVec, 1, 0, max, sum);
867 metrics[next]->setName(
"edge cuts");
868 metrics[next]->setMetricValue(
"global maximum", max);
869 metrics[next]->setMetricValue(
"global sum", sum);
874 scalar_t *wgt = sumBuf + nparts;
876 for (
int edim=0; edim < ewgtDim; edim++){
877 ArrayView<scalar_t> fromVec(wgt, nparts);
878 getStridedStats<scalar_t>(fromVec, 1, 0, max, sum);
880 std::ostringstream oss;
881 oss <<
"weight " << edim;
883 metrics[next]->setName(oss.str());
884 metrics[next]->setMetricValue(
"global maximum", max);
885 metrics[next]->setMetricValue(
"global sum", sum);
899 template <
typename scalar_t,
typename part_t>
902 os <<
"Graph Metrics: (" << numParts <<
" parts)";
904 if (targetNumParts != numParts) {
905 os <<
"Target number of parts is: " << targetNumParts << std::endl;
912 template <
typename scalar_t,
typename part_t>
915 printGraphMetricsHeader<scalar_t, part_t>(os, targetNumParts, numParts);
916 for (
int i=0; i < infoList.size(); i++) {
918 infoList[i]->printLine(os);
926 template <
typename scalar_t,
typename part_t>
929 printGraphMetricsHeader<scalar_t, part_t>(os, targetNumParts, numParts);
930 metricValue->printLine(os);