13 #ifndef ZOLTAN2_GRAPHICMETRICVALUESUTILITY_HPP
14 #define ZOLTAN2_GRAPHICMETRICVALUESUTILITY_HPP
19 #include <zoltan_dd.h>
54 template <
typename Adapter,
56 MachineRepresentation<
typename Adapter::scalar_t,
59 const RCP<const Environment> &env,
60 const RCP<
const Comm<int> > &comm,
62 const ArrayView<const typename Adapter::part_t> &parts,
65 ArrayRCP<typename Adapter::scalar_t> &globalSums,
66 bool bMessages =
true,
67 const RCP <const MachineRep> machine = Teuchos::null) {
74 bool bHops = (machine != Teuchos::null);
83 typedef typename Adapter::offset_t offset_t;
84 typedef typename Adapter::scalar_t scalar_t;
91 lno_t localNumEdges = graph->getLocalNumEdges();
93 ArrayView<const gno_t> Ids;
94 ArrayView<t_input_t> v_wghts;
95 graph->getVertexList(Ids, v_wghts);
100 ArrayView<const gno_t> edgeIds;
101 ArrayView<const offset_t> offsets;
102 ArrayView<t_input_t> e_wgts;
103 graph->getEdgeList(edgeIds, offsets, e_wgts);
106 std::vector <scalar_t> edge_weights;
107 int numWeightPerEdge = graph->getNumWeightsPerEdge();
109 int numMetrics = bHops ?
113 if (numWeightPerEdge) numMetrics += bHops ?
114 numWeightPerEdge * 2 :
118 auto next = metrics.size();
119 for (
auto n = 0; n < numMetrics; ++n) {
120 addNewMetric<gm_t, scalar_t>(env, metrics);
123 std::vector <part_t> e_parts (localNumEdges);
125 std::vector<part_t> local_parts;
127 #ifdef HAVE_ZOLTAN2_MPI
128 if (comm->getSize() > 1) {
129 const bool bUseLocalIDs =
false;
132 directory_t directory(comm, bUseLocalIDs, debug_level);
134 if (localNumVertices)
135 directory.update(localNumVertices, &Ids[0], NULL, &parts[0],
138 directory.update(localNumVertices, NULL, NULL, NULL,
142 directory.find(localNumEdges, &edgeIds[0], NULL, &e_parts[0],
145 directory.find(localNumEdges, NULL, NULL, NULL,
150 std::map<gno_t,lno_t> global_id_to_local_index;
156 for (lno_t i = 0; i < localNumVertices; ++i){
159 global_id_to_local_index[Ids[i]] = i;
162 for (lno_t i = 0; i < localNumEdges; ++i){
163 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");
174 const bool bUseLocalIDs =
false;
182 part_info() : sum_weights(0) {
189 const part_info & operator+=(
const part_info & src) {
190 sum_weights += src.sum_weights;
196 bool operator>(
const part_info & src) {
208 return (target_part > src.target_part);
219 bool operator==(
const part_info & src) {
222 return (target_part == src.target_part);
226 scalar_t sum_weights;
230 std::vector <lno_t> part_begins(numParts, -1);
231 std::vector <lno_t> part_nexts(localNumVertices, -1);
235 for (lno_t i = 0; i < localNumVertices; ++i){
236 part_t ap = parts[i];
237 part_nexts[i] = part_begins[ap];
241 for (
int weight_index = -1;
242 weight_index < numWeightPerEdge; ++weight_index) {
244 std::vector<part_t> part_data(numParts);
245 std::vector<std::vector<part_info>> user_data(numParts);
246 int count_total_entries = 0;
248 std::vector <part_t> part_neighbors(numParts);
249 std::vector <scalar_t> part_neighbor_weights(numParts, 0);
250 std::vector <scalar_t> part_neighbor_weights_ordered(numParts);
253 for (part_t i = 0; i < numParts; ++i) {
254 part_t num_neighbor_parts = 0;
255 lno_t v = part_begins[i];
259 for (offset_t j = offsets[v]; j < offsets[v+1]; ++j){
262 part_t ep = e_parts[j];
268 if (weight_index > -1){
269 ew = e_wgts[weight_index][j];
272 if (part_neighbor_weights[ep] < 0.00001){
273 part_neighbors[num_neighbor_parts++] = ep;
275 part_neighbor_weights[ep] += ew;
281 for (lno_t j = 0; j < num_neighbor_parts; ++j) {
282 part_t neighbor_part = part_neighbors[j];
283 part_neighbor_weights_ordered[j] =
284 part_neighbor_weights[neighbor_part];
285 part_neighbor_weights[neighbor_part] = 0;
288 if (num_neighbor_parts > 0) {
294 part_data[count_total_entries] = i;
295 std::vector<part_info> & add_user_data =
296 user_data[count_total_entries];
297 ++count_total_entries;
299 add_user_data.resize(num_neighbor_parts);
301 for(
int n = 0; n < num_neighbor_parts; ++n) {
302 part_info & add_data = add_user_data[n];
303 add_data.target_part = part_neighbors[n];
304 add_data.sum_weights = part_neighbor_weights_ordered[n];
309 scalar_t max_edge_cut = 0;
310 scalar_t total_edge_cut = 0;
311 part_t max_message = 0;
312 part_t total_message = 0;
314 part_t total_hop_count = 0;
315 scalar_t total_weighted_hop_count = 0;
316 part_t max_hop_count = 0;
317 scalar_t max_weighted_hop_count = 0;
322 if(local_parts.size() == 0) {
323 local_parts.resize(numParts);
324 for(
size_t n = 0; n < local_parts.size(); ++n) {
329 std::vector<std::vector<part_info>> find_user_data;
338 (
dynamic_cast<const Teuchos::SerialComm<int>*
>(&(*comm)) != NULL);
343 directory_t directory(comm, bUseLocalIDs, debug_level);
345 if(count_total_entries) {
347 directory.update(count_total_entries, &part_data[0],
349 NULL, directory_t::Update_Mode::AggregateAdd);
352 directory.update(count_total_entries, NULL, NULL, NULL,
353 NULL, directory_t::Update_Mode::AggregateAdd);
357 directory.get_locally_managed_gids(local_parts);
360 find_user_data.resize(local_parts.size());
363 directory.find(local_parts.size(), &local_parts[0], NULL,
364 &find_user_data[0], NULL, NULL,
false);
367 find_user_data = user_data;
370 for(
size_t n = 0; n < local_parts.size(); ++n) {
371 scalar_t part_edge_cut = 0;
372 part_t part_messages = 0;
373 const std::vector<part_info> & data = find_user_data[n];
374 for(
size_t q = 0; q < data.size(); ++q) {
375 const part_t local_part = local_parts[n];
376 const part_info & info = data[q];
377 if (info.target_part != local_part) {
378 part_edge_cut += info.sum_weights;
382 typename MachineRep::machine_pcoord_t hop_count = 0;
383 machine->getHopCount(local_part, info.target_part, hop_count);
385 part_t hop_counts = hop_count;
386 scalar_t weighted_hop_counts = hop_count * info.sum_weights;
388 total_hop_count += hop_counts;
389 total_weighted_hop_count += weighted_hop_counts;
391 if (hop_counts > max_hop_count ){
392 max_hop_count = hop_counts;
394 if (weighted_hop_counts > max_weighted_hop_count ){
395 max_weighted_hop_count = weighted_hop_counts;
401 if(part_edge_cut > max_edge_cut) {
402 max_edge_cut = part_edge_cut;
404 total_edge_cut += part_edge_cut;
406 if (part_messages > max_message){
407 max_message = part_messages;
409 total_message += part_messages;
413 scalar_t g_max_edge_cut = 0;
414 scalar_t g_total_edge_cut = 0;
415 part_t g_max_message = 0;
416 part_t g_total_message = 0;
418 part_t g_total_hop_count = 0;
419 scalar_t g_total_weighted_hop_count = 0;
420 part_t g_max_hop_count = 0;
421 scalar_t g_max_weighted_hop_count = 0;
424 Teuchos::reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_MAX, 1,
425 &max_edge_cut, &g_max_edge_cut);
426 Teuchos::reduceAll<int, part_t>(*comm, Teuchos::REDUCE_MAX, 1,
427 &max_message, &g_max_message);
429 Teuchos::reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_SUM, 1,
430 &total_edge_cut, &g_total_edge_cut);
431 Teuchos::reduceAll<int, part_t>(*comm, Teuchos::REDUCE_SUM, 1,
432 &total_message, &g_total_message);
435 Teuchos::reduceAll<int, part_t>(*comm, Teuchos::REDUCE_MAX, 1,
436 &max_hop_count, &g_max_hop_count);
437 Teuchos::reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_MAX, 1,
438 &max_weighted_hop_count,
439 &g_max_weighted_hop_count);
441 Teuchos::reduceAll<int, part_t>(*comm, Teuchos::REDUCE_SUM, 1,
442 &total_hop_count, &g_total_hop_count);
443 Teuchos::reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_SUM, 1,
444 &total_weighted_hop_count,
445 &g_total_weighted_hop_count);
450 if (weight_index == -1){
451 metrics[next]->setName(
"edge cuts");
454 std::ostringstream oss;
455 oss <<
"weight " << weight_index;
456 metrics[next]->setName( oss.str());
458 metrics[next]->setMetricValue(
"global maximum", g_max_edge_cut);
459 metrics[next]->setMetricValue(
"global sum", g_total_edge_cut);
461 if (weight_index == -1){
462 metrics[next]->setName(
"message");
463 metrics[next]->setMetricValue(
"global maximum", g_max_message);
464 metrics[next]->setMetricValue(
"global sum", g_total_message);
469 if (weight_index == -1){
470 metrics[next]->setName(
"hops (No Weight)");
471 metrics[next]->setMetricValue(
"global maximum", g_max_hop_count);
472 metrics[next]->setMetricValue(
"global sum", g_total_hop_count);
476 std::ostringstream oss;
477 oss <<
"weighted hops" << weight_index;
478 metrics[next]->setName( oss.str());
480 setMetricValue(
"global maximum", g_max_weighted_hop_count);
482 setMetricValue(
"global sum", g_total_weighted_hop_count);
495 template <
typename scalar_t,
typename part_t>
500 os <<
"Graph Metrics: (" << numParts <<
" parts)";
502 if (targetNumParts != numParts) {
503 os <<
"Target number of parts is: " << targetNumParts << std::endl;
510 template <
typename scalar_t,
typename part_t>
517 printGraphMetricsHeader<scalar_t, part_t>(os, targetNumParts, numParts);
518 for (
int i=0; i < infoList.size(); i++) {
520 infoList[i]->printLine(os);
528 template <
typename scalar_t,
typename part_t>
534 printGraphMetricsHeader<scalar_t, part_t>(os, targetNumParts, numParts);
535 metricValue->printLine(os);
Time an algorithm (or other entity) as a whole.
map_t::global_ordinal_type gno_t
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.
sub-steps, each method's entry and exit
SparseMatrixAdapter_t::part_t part_t
void globalWeightedByPart(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, bool bMessages=true, const RCP< const MachineRep > machine=Teuchos::null)
Given the local partitioning, compute the global weighted cuts in each part.
map_t::local_ordinal_type lno_t
void printGraphMetrics(std::ostream &os, part_t targetNumParts, part_t numParts, const ArrayView< RCP< BaseClassMetrics< scalar_t > > > &infoList)
Print out list of graph metrics.
size_t getLocalNumVertices() const
Returns the number vertices on this process.
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.
#define Z2_THROW_OUTSIDE_ERROR(env)
Throw an error returned from outside the Zoltan2 library.
#define METRICS_UNSET_STRING