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::scalar_t t_scalar_t;
81 typedef typename Adapter::node_t t_node_t;
87 t_lno_t localNumEdges = graph->getLocalNumEdges();
89 ArrayView<const t_gno_t> Ids;
90 ArrayView<t_input_t> v_wghts;
91 graph->getVertexList(Ids, v_wghts);
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);
102 std::vector <t_scalar_t> edge_weights;
103 int numWeightPerEdge = graph->getNumWeightsPerEdge();
106 if (numWeightPerEdge) numMetrics += numWeightPerEdge * 2;
109 auto next = metrics.size();
110 for (
auto n = 0; n < numMetrics; ++n) {
111 addNewMetric<gm_t, t_scalar_t>(env, metrics);
114 std::vector <part_t> e_parts (localNumEdges);
115 #ifdef HAVE_ZOLTAN2_MPI 116 if (comm->getSize() > 1)
118 Zoltan_DD_Struct *dd = NULL;
120 MPI_Comm mpicomm = Teuchos::getRawMpiComm(*comm);
124 Zoltan_DD_Create(&dd, mpicomm,
126 sizeof(part_t), localNumVertices, debug_level);
128 ZOLTAN_ID_PTR ddnotneeded = NULL;
131 (localNumVertices ? (ZOLTAN_ID_PTR) Ids.getRawPtr() : NULL),
133 (localNumVertices ? (
char *) &(parts[0]) : NULL),
135 int(localNumVertices));
139 (localNumEdges ? (ZOLTAN_ID_PTR) edgeIds.getRawPtr() : NULL),
141 (localNumEdges ? (
char *)&(e_parts[0]) : NULL),
146 Zoltan_DD_Destroy(&dd);
151 std::map<t_gno_t,t_lno_t> global_id_to_local_index;
156 for (t_lno_t i = 0; i < localNumVertices; ++i){
159 global_id_to_local_index[Ids[i]] = i;
162 for (t_lno_t i = 0; i < localNumEdges; ++i){
163 t_gno_t ei = edgeIds[i];
165 part_t p = parts[global_id_to_local_index[ei]];
170 RCP<const Teuchos::Comm<int> > tcomm = comm;
172 env->timerStart(
MACRO_TIMERS,
"Communication Graph Create");
175 std::vector <t_lno_t> part_begins(numParts, -1);
176 std::vector <t_lno_t> part_nexts(localNumVertices, -1);
180 for (t_lno_t i = 0; i < localNumVertices; ++i){
181 part_t ap = parts[i];
182 part_nexts[i] = part_begins[ap];
187 for (
int weight_index = -1; weight_index < numWeightPerEdge ; ++weight_index){
192 typedef t_lno_t local_part_type;
193 typedef t_gno_t global_part_type;
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));
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));
202 std::vector <global_part_type> part_neighbors (numParts);
204 std::vector <t_scalar_t> part_neighbor_weights(numParts, 0);
205 std::vector <t_scalar_t> part_neighbor_weights_ordered(numParts);
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];
214 for (t_lno_t j = offsets[v]; j < offsets[v+1]; ++j){
216 part_t ep = e_parts[j];
219 if (weight_index > -1){
220 ew = e_wgts[weight_index][j];
223 if (part_neighbor_weights[ep] < 0.00001){
224 part_neighbors[num_neighbor_parts++] = ep;
226 part_neighbor_weights[ep] += ew;
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;
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);
245 tMatrix->fillComplete ();
246 local_part_type num_local_parts = map->getNodeNumElements();
248 Array<global_part_type> Indices;
249 Array<t_scalar_t> Values;
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;
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;
261 for (local_part_type i=0; i < num_local_parts; i++) {
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);
269 t_scalar_t part_edge_cut = 0;
270 global_part_type part_messages = 0;
272 for (
size_t j=0; j < NumEntries; j++){
273 if (Indices[j] != globalRow){
274 part_edge_cut += Values[j];
277 typename MachineRep::machine_pcoord_t hop_count = 0;
278 machine->getHopCount(globalRow, Indices[j], hop_count);
280 global_part_type hop_counts = hop_count;
281 t_scalar_t weighted_hop_counts = hop_count * Values[j];
283 total_hop_count += hop_counts;
284 total_weighted_hop_count += weighted_hop_counts;
286 if (hop_counts > max_hop_count ){
287 max_hop_count = hop_counts;
289 if (weighted_hop_counts > max_weighted_hop_count ){
290 max_weighted_hop_count = weighted_hop_counts;
294 if (part_edge_cut > max_edge_cut){
295 max_edge_cut = part_edge_cut;
297 total_edge_cut += part_edge_cut;
299 if (part_messages > max_message){
300 max_message = part_messages;
302 total_message += part_messages;
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;
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;
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);
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);
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);
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);
335 if (weight_index == -1){
336 metrics[next]->setName(
"md edge cuts");
339 std::ostringstream oss;
340 oss <<
"md weighted edge cuts" << weight_index;
341 metrics[next]->setName( oss.str());
344 metrics[next]->setMetricValue(
"global maximum", g_max_edge_cut);
345 metrics[next]->setMetricValue(
"global sum", g_total_edge_cut);
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);
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);
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);
372 env->timerStop(
MACRO_TIMERS,
"Communication Graph Create");
374 env->debug(
DETAILED_STATUS,
"Exiting globalWeightedCutsMessagesHopsByPart");
378 template <
typename Adapter>
380 const RCP<const Environment> &env,
381 const RCP<
const Comm<int> > &comm,
383 const ArrayView<const typename Adapter::part_t> &parts,
386 ArrayRCP<typename Adapter::scalar_t> &globalSums)
388 env->debug(
DETAILED_STATUS,
"Entering globalWeightedCutsMessagesByPart");
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;
396 typedef typename Adapter::node_t t_node_t;
402 t_lno_t localNumEdges = graph->getLocalNumEdges();
404 ArrayView<const t_gno_t> Ids;
405 ArrayView<t_input_t> v_wghts;
406 graph->getVertexList(Ids, v_wghts);
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);
417 std::vector <t_scalar_t> edge_weights;
418 int numWeightPerEdge = graph->getNumWeightsPerEdge();
421 if (numWeightPerEdge) numMetrics += numWeightPerEdge;
424 auto next = metrics.size();
425 for (
auto n = 0; n < numMetrics; ++n) {
426 addNewMetric<gm_t, t_scalar_t>(env, metrics);
429 std::vector <part_t> e_parts (localNumEdges);
430 #ifdef HAVE_ZOLTAN2_MPI 431 if (comm->getSize() > 1)
433 Zoltan_DD_Struct *dd = NULL;
435 MPI_Comm mpicomm = Teuchos::getRawMpiComm(*comm);
439 Zoltan_DD_Create(&dd, mpicomm,
441 sizeof(part_t), localNumVertices, debug_level);
443 ZOLTAN_ID_PTR ddnotneeded = NULL;
446 (localNumVertices ? (ZOLTAN_ID_PTR) Ids.getRawPtr() : NULL),
448 (localNumVertices ? (
char *) &(parts[0]) : NULL),
450 int(localNumVertices));
454 (localNumEdges ? (ZOLTAN_ID_PTR) edgeIds.getRawPtr() : NULL),
456 (localNumEdges ? (
char *)&(e_parts[0]) : NULL),
461 Zoltan_DD_Destroy(&dd);
466 std::map<t_gno_t,t_lno_t> global_id_to_local_index;
471 for (t_lno_t i = 0; i < localNumVertices; ++i){
474 global_id_to_local_index[Ids[i]] = i;
477 for (t_lno_t i = 0; i < localNumEdges; ++i){
478 t_gno_t ei = edgeIds[i];
480 part_t p = parts[global_id_to_local_index[ei]];
485 RCP<const Teuchos::Comm<int> > tcomm = comm;
487 env->timerStart(
MACRO_TIMERS,
"Communication Graph Create");
490 std::vector <t_lno_t> part_begins(numParts, -1);
491 std::vector <t_lno_t> part_nexts(localNumVertices, -1);
495 for (t_lno_t i = 0; i < localNumVertices; ++i){
496 part_t ap = parts[i];
497 part_nexts[i] = part_begins[ap];
501 for (
int weight_index = -1; weight_index < numWeightPerEdge ; ++weight_index){
506 typedef t_lno_t local_part_type;
507 typedef t_gno_t global_part_type;
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));
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));
516 std::vector <global_part_type> part_neighbors (numParts);
518 std::vector <t_scalar_t> part_neighbor_weights(numParts, 0);
519 std::vector <t_scalar_t> part_neighbor_weights_ordered(numParts);
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];
528 for (t_lno_t j = offsets[v]; j < offsets[v+1]; ++j){
530 part_t ep = e_parts[j];
533 if (weight_index > -1){
534 ew = e_wgts[weight_index][j];
537 if (part_neighbor_weights[ep] < 0.00001){
538 part_neighbors[num_neighbor_parts++] = ep;
540 part_neighbor_weights[ep] += ew;
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;
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);
559 tMatrix->fillComplete ();
560 local_part_type num_local_parts = map->getNodeNumElements();
562 Array<global_part_type> Indices;
563 Array<t_scalar_t> Values;
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;
570 for (local_part_type i=0; i < num_local_parts; i++) {
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);
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];
586 if (part_edge_cut > max_edge_cut){
587 max_edge_cut = part_edge_cut;
589 total_edge_cut += part_edge_cut;
591 if (part_messages > max_message){
592 max_message = part_messages;
594 total_message += part_messages;
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;
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);
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);
612 if (weight_index == -1){
613 metrics[next]->setName(
"md edge cuts");
616 std::ostringstream oss;
617 oss <<
"md weight " << weight_index;
618 metrics[next]->setName( oss.str());
620 metrics[next]->setMetricValue(
"global maximum", g_max_edge_cut);
621 metrics[next]->setMetricValue(
"global sum", g_total_edge_cut);
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);
631 env->timerStop(
MACRO_TIMERS,
"Communication Graph Create");
633 env->debug(
DETAILED_STATUS,
"Exiting globalWeightedCutsMessagesByPart");
664 template <
typename Adapter>
666 const RCP<const Environment> &env,
667 const RCP<
const Comm<int> > &comm,
669 const ArrayView<const typename Adapter::part_t> &part,
672 ArrayRCP<typename Adapter::scalar_t> &globalSums)
680 int ewgtDim = graph->getNumWeightsPerEdge();
683 if (ewgtDim) numMetrics += ewgtDim;
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;
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;
696 const GST
INVALID = Teuchos::OrdinalTraits<GST>::invalid ();
700 auto next = metrics.size();
702 for (
auto n = 0; n < numMetrics; ++n) {
703 addNewMetric<gm_t, scalar_t>(env, metrics);
710 lno_t localNumObj = part.size();
711 part_t localNum[2], globalNum[2];
712 localNum[0] =
static_cast<part_t
>(ewgtDim);
715 for (lno_t i=0; i < localNumObj; i++)
716 if (part[i] > localNum[1]) localNum[1] = part[i];
719 reduceAll<int, part_t>(*comm, Teuchos::REDUCE_MAX, 2,
720 localNum, globalNum);
724 env->globalBugAssertion(__FILE__,__LINE__,
725 "inconsistent number of edge weights",
728 part_t nparts = globalNum[1] + 1;
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);
738 scalar_t *localBuf =
new scalar_t [globalSumSize];
739 env->localMemoryAssertion(__FILE__,__LINE__,globalSumSize,localBuf);
740 memset(localBuf, 0,
sizeof(scalar_t) * globalSumSize);
742 scalar_t *cut = localBuf;
744 ArrayView<const gno_t> Ids;
745 ArrayView<input_t> vwgts;
747 graph->getVertexList(Ids, vwgts);
749 ArrayView<const gno_t> edgeIds;
750 ArrayView<const lno_t> offsets;
751 ArrayView<input_t> wgts;
753 graph->getEdgeList(edgeIds, offsets, wgts);
758 RCP<const map_type> vertexMapG;
761 gno_t min = std::numeric_limits<gno_t>::max();
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;
770 Teuchos::reduceAll<int, gno_t>(*comm,Teuchos::REDUCE_MIN,1,&min,&gmin);
773 vertexMapG = rcp(
new map_type(INVALID, Ids, gmin, comm));
781 RCP<sparse_matrix_type> adjsMatrix;
784 adjsMatrix = rcp (
new sparse_matrix_type (vertexMapG, 0));
786 Array<part_t> justOneA(maxcols, 1);
788 for (lno_t localElement=0; localElement<localNumObj; ++localElement){
790 size_t ncols = offsets[localElement+1] - offsets[localElement];
791 adjsMatrix->insertGlobalValues(Ids[localElement],
792 edgeIds(offsets[localElement], ncols),
797 adjsMatrix->fillComplete ();
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]);
806 adjsMatrix->rightScale(*scaleVec);
807 Array<gno_t> Indices;
808 Array<part_t> Values;
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);
817 for (
size_t j=0; j < NumEntries; j++)
818 if (part[i] != Values[j])
822 if (numMetrics > 1) {
824 scalar_t *wgt = localBuf + nparts;
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);
836 for (
size_t j=0; j < NumEntries; j++)
837 if (part[i] != Values[j])
838 wgt[part[i]] += wgts[edim][offsets[i] + j];
848 reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_SUM, globalSumSize,
859 scalar_t max=0, sum=0;
861 ArrayView<scalar_t> cutVec(cut, nparts);
862 getStridedStats<scalar_t>(cutVec, 1, 0, max, sum);
864 metrics[next]->setName(
"edge cuts");
865 metrics[next]->setMetricValue(
"global maximum", max);
866 metrics[next]->setMetricValue(
"global sum", sum);
871 scalar_t *wgt = sumBuf + nparts;
873 for (
int edim=0; edim < ewgtDim; edim++){
874 ArrayView<scalar_t> fromVec(wgt, nparts);
875 getStridedStats<scalar_t>(fromVec, 1, 0, max, sum);
877 std::ostringstream oss;
878 oss <<
"weight " << edim;
880 metrics[next]->setName(oss.str());
881 metrics[next]->setMetricValue(
"global maximum", max);
882 metrics[next]->setMetricValue(
"global sum", sum);
896 template <
typename scalar_t,
typename part_t>
899 os <<
"Graph Metrics: (" << numParts <<
" parts)";
901 if (targetNumParts != numParts) {
902 os <<
"Target number of parts is: " << targetNumParts << std::endl;
909 template <
typename scalar_t,
typename part_t>
912 printGraphMetricsHeader<scalar_t, part_t>(os, targetNumParts, numParts);
913 for (
int i=0; i < infoList.size(); i++) {
915 infoList[i]->printLine(os);
923 template <
typename scalar_t,
typename part_t>
926 printGraphMetricsHeader<scalar_t, part_t>(os, targetNumParts, numParts);
927 metricValue->printLine(os);
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tcrsMatrix_t
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'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
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'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