49 #ifndef ZOLTAN2_GRAPHICMETRICVALUESUTILITY_HPP
50 #define ZOLTAN2_GRAPHICMETRICVALUESUTILITY_HPP
55 #include <zoltan_dd.h>
90 template <
typename Adapter,
92 MachineRepresentation<
typename Adapter::scalar_t,
95 const RCP<const Environment> &env,
96 const RCP<
const Comm<int> > &comm,
98 const ArrayView<const typename Adapter::part_t> &parts,
101 ArrayRCP<typename Adapter::scalar_t> &globalSums,
102 bool bMessages =
true,
103 const RCP <const MachineRep> machine = Teuchos::null) {
110 bool bHops = (machine != Teuchos::null);
117 typedef typename Adapter::lno_t lno_t;
118 typedef typename Adapter::gno_t gno_t;
119 typedef typename Adapter::offset_t offset_t;
120 typedef typename Adapter::scalar_t scalar_t;
127 lno_t localNumEdges = graph->getLocalNumEdges();
129 ArrayView<const gno_t> Ids;
130 ArrayView<t_input_t> v_wghts;
131 graph->getVertexList(Ids, v_wghts);
136 ArrayView<const gno_t> edgeIds;
137 ArrayView<const offset_t> offsets;
138 ArrayView<t_input_t> e_wgts;
139 graph->getEdgeList(edgeIds, offsets, e_wgts);
142 std::vector <scalar_t> edge_weights;
143 int numWeightPerEdge = graph->getNumWeightsPerEdge();
145 int numMetrics = bHops ?
149 if (numWeightPerEdge) numMetrics += bHops ?
150 numWeightPerEdge * 2 :
154 auto next = metrics.size();
155 for (
auto n = 0; n < numMetrics; ++n) {
156 addNewMetric<gm_t, scalar_t>(env, metrics);
159 std::vector <part_t> e_parts (localNumEdges);
161 std::vector<part_t> local_parts;
163 #ifdef HAVE_ZOLTAN2_MPI
164 if (comm->getSize() > 1) {
165 const bool bUseLocalIDs =
false;
168 directory_t directory(comm, bUseLocalIDs, debug_level);
170 if (localNumVertices)
171 directory.update(localNumVertices, &Ids[0], NULL, &parts[0],
172 NULL, directory_t::Update_Mode::Replace);
174 directory.update(localNumVertices, NULL, NULL, NULL,
175 NULL, directory_t::Update_Mode::Replace);
178 directory.find(localNumEdges, &edgeIds[0], NULL, &e_parts[0],
181 directory.find(localNumEdges, NULL, NULL, NULL,
186 std::map<gno_t,lno_t> global_id_to_local_index;
192 for (lno_t i = 0; i < localNumVertices; ++i){
195 global_id_to_local_index[Ids[i]] = i;
198 for (lno_t i = 0; i < localNumEdges; ++i){
199 gno_t ei = edgeIds[i];
201 part_t p = parts[global_id_to_local_index[ei]];
206 RCP<const Teuchos::Comm<int> > tcomm = comm;
208 env->timerStart(
MACRO_TIMERS,
"Communication Graph Create");
210 const bool bUseLocalIDs =
false;
218 part_info() : sum_weights(0) {
225 const part_info & operator+=(
const part_info & src) {
226 sum_weights += src.sum_weights;
232 bool operator>(
const part_info & src) {
244 return (target_part > src.target_part);
255 bool operator==(
const part_info & src) {
258 return (target_part == src.target_part);
262 scalar_t sum_weights;
266 std::vector <lno_t> part_begins(numParts, -1);
267 std::vector <lno_t> part_nexts(localNumVertices, -1);
271 for (lno_t i = 0; i < localNumVertices; ++i){
272 part_t ap = parts[i];
273 part_nexts[i] = part_begins[ap];
277 for (
int weight_index = -1;
278 weight_index < numWeightPerEdge; ++weight_index) {
280 std::vector<part_t> part_data(numParts);
281 std::vector<std::vector<part_info>> user_data(numParts);
282 int count_total_entries = 0;
284 std::vector <part_t> part_neighbors(numParts);
285 std::vector <scalar_t> part_neighbor_weights(numParts, 0);
286 std::vector <scalar_t> part_neighbor_weights_ordered(numParts);
289 for (part_t i = 0; i < numParts; ++i) {
290 part_t num_neighbor_parts = 0;
291 lno_t v = part_begins[i];
295 for (offset_t j = offsets[v]; j < offsets[v+1]; ++j){
298 part_t ep = e_parts[j];
304 if (weight_index > -1){
305 ew = e_wgts[weight_index][j];
308 if (part_neighbor_weights[ep] < 0.00001){
309 part_neighbors[num_neighbor_parts++] = ep;
311 part_neighbor_weights[ep] += ew;
317 for (lno_t j = 0; j < num_neighbor_parts; ++j) {
318 part_t neighbor_part = part_neighbors[j];
319 part_neighbor_weights_ordered[j] =
320 part_neighbor_weights[neighbor_part];
321 part_neighbor_weights[neighbor_part] = 0;
324 if (num_neighbor_parts > 0) {
330 part_data[count_total_entries] = i;
331 std::vector<part_info> & add_user_data =
332 user_data[count_total_entries];
333 ++count_total_entries;
335 add_user_data.resize(num_neighbor_parts);
337 for(
int n = 0; n < num_neighbor_parts; ++n) {
338 part_info & add_data = add_user_data[n];
339 add_data.target_part = part_neighbors[n];
340 add_data.sum_weights = part_neighbor_weights_ordered[n];
345 scalar_t max_edge_cut = 0;
346 scalar_t total_edge_cut = 0;
347 part_t max_message = 0;
348 part_t total_message = 0;
350 part_t total_hop_count = 0;
351 scalar_t total_weighted_hop_count = 0;
352 part_t max_hop_count = 0;
353 scalar_t max_weighted_hop_count = 0;
358 if(local_parts.size() == 0) {
359 local_parts.resize(numParts);
360 for(
size_t n = 0; n < local_parts.size(); ++n) {
365 std::vector<std::vector<part_info>> find_user_data;
374 (
dynamic_cast<const Teuchos::SerialComm<int>*
>(&(*comm)) != NULL);
379 directory_t directory(comm, bUseLocalIDs, debug_level);
381 if(count_total_entries) {
383 directory.update(count_total_entries, &part_data[0],
385 NULL, directory_t::Update_Mode::AggregateAdd);
388 directory.update(count_total_entries, NULL, NULL, NULL,
389 NULL, directory_t::Update_Mode::AggregateAdd);
393 directory.get_locally_managed_gids(local_parts);
396 find_user_data.resize(local_parts.size());
399 directory.find(local_parts.size(), &local_parts[0], NULL,
400 &find_user_data[0], NULL, NULL,
false);
403 find_user_data = user_data;
406 for(
size_t n = 0; n < local_parts.size(); ++n) {
407 scalar_t part_edge_cut = 0;
408 part_t part_messages = 0;
409 const std::vector<part_info> & data = find_user_data[n];
410 for(
size_t q = 0; q < data.size(); ++q) {
411 const part_t local_part = local_parts[n];
412 const part_info & info = data[q];
413 if (info.target_part != local_part) {
414 part_edge_cut += info.sum_weights;
418 typename MachineRep::machine_pcoord_t hop_count = 0;
419 machine->getHopCount(local_part, info.target_part, hop_count);
421 part_t hop_counts = hop_count;
422 scalar_t weighted_hop_counts = hop_count * info.sum_weights;
424 total_hop_count += hop_counts;
425 total_weighted_hop_count += weighted_hop_counts;
427 if (hop_counts > max_hop_count ){
428 max_hop_count = hop_counts;
430 if (weighted_hop_counts > max_weighted_hop_count ){
431 max_weighted_hop_count = weighted_hop_counts;
437 if(part_edge_cut > max_edge_cut) {
438 max_edge_cut = part_edge_cut;
440 total_edge_cut += part_edge_cut;
442 if (part_messages > max_message){
443 max_message = part_messages;
445 total_message += part_messages;
449 scalar_t g_max_edge_cut = 0;
450 scalar_t g_total_edge_cut = 0;
451 part_t g_max_message = 0;
452 part_t g_total_message = 0;
454 part_t g_total_hop_count = 0;
455 scalar_t g_total_weighted_hop_count = 0;
456 part_t g_max_hop_count = 0;
457 scalar_t g_max_weighted_hop_count = 0;
460 Teuchos::reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_MAX, 1,
461 &max_edge_cut, &g_max_edge_cut);
462 Teuchos::reduceAll<int, part_t>(*comm, Teuchos::REDUCE_MAX, 1,
463 &max_message, &g_max_message);
465 Teuchos::reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_SUM, 1,
466 &total_edge_cut, &g_total_edge_cut);
467 Teuchos::reduceAll<int, part_t>(*comm, Teuchos::REDUCE_SUM, 1,
468 &total_message, &g_total_message);
471 Teuchos::reduceAll<int, part_t>(*comm, Teuchos::REDUCE_MAX, 1,
472 &max_hop_count, &g_max_hop_count);
473 Teuchos::reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_MAX, 1,
474 &max_weighted_hop_count,
475 &g_max_weighted_hop_count);
477 Teuchos::reduceAll<int, part_t>(*comm, Teuchos::REDUCE_SUM, 1,
478 &total_hop_count, &g_total_hop_count);
479 Teuchos::reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_SUM, 1,
480 &total_weighted_hop_count,
481 &g_total_weighted_hop_count);
486 if (weight_index == -1){
487 metrics[next]->setName(
"edge cuts");
490 std::ostringstream oss;
491 oss <<
"weight " << weight_index;
492 metrics[next]->setName( oss.str());
494 metrics[next]->setMetricValue(
"global maximum", g_max_edge_cut);
495 metrics[next]->setMetricValue(
"global sum", g_total_edge_cut);
497 if (weight_index == -1){
498 metrics[next]->setName(
"message");
499 metrics[next]->setMetricValue(
"global maximum", g_max_message);
500 metrics[next]->setMetricValue(
"global sum", g_total_message);
505 if (weight_index == -1){
506 metrics[next]->setName(
"hops (No Weight)");
507 metrics[next]->setMetricValue(
"global maximum", g_max_hop_count);
508 metrics[next]->setMetricValue(
"global sum", g_total_hop_count);
512 std::ostringstream oss;
513 oss <<
"weighted hops" << weight_index;
514 metrics[next]->setName( oss.str());
516 setMetricValue(
"global maximum", g_max_weighted_hop_count);
518 setMetricValue(
"global sum", g_total_weighted_hop_count);
531 template <
typename scalar_t,
typename part_t>
536 os <<
"Graph Metrics: (" << numParts <<
" parts)";
538 if (targetNumParts != numParts) {
539 os <<
"Target number of parts is: " << targetNumParts << std::endl;
546 template <
typename scalar_t,
typename part_t>
553 printGraphMetricsHeader<scalar_t, part_t>(os, targetNumParts, numParts);
554 for (
int i=0; i < infoList.size(); i++) {
556 infoList[i]->printLine(os);
564 template <
typename scalar_t,
typename part_t>
570 printGraphMetricsHeader<scalar_t, part_t>(os, targetNumParts, numParts);
571 metricValue->printLine(os);
Time an algorithm (or other entity) as a whole.
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.
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