ViennaCL - The Vienna Computing Library  1.5.0
viennacl/misc/cuthill_mckee.hpp
Go to the documentation of this file.
00001 #ifndef VIENNACL_MISC_CUTHILL_MCKEE_HPP
00002 #define VIENNACL_MISC_CUTHILL_MCKEE_HPP
00003 
00004 /* =========================================================================
00005    Copyright (c) 2010-2013, Institute for Microelectronics,
00006                             Institute for Analysis and Scientific Computing,
00007                             TU Wien.
00008    Portions of this software are copyright by UChicago Argonne, LLC.
00009 
00010                             -----------------
00011                   ViennaCL - The Vienna Computing Library
00012                             -----------------
00013 
00014    Project Head:    Karl Rupp                   rupp@iue.tuwien.ac.at
00015 
00016    (A list of authors and contributors can be found in the PDF manual)
00017 
00018    License:         MIT (X11), see file LICENSE in the base directory
00019 ============================================================================= */
00020 
00021 
00028 #include <iostream>
00029 #include <iterator>
00030 #include <fstream>
00031 #include <string>
00032 #include <algorithm>
00033 #include <map>
00034 #include <vector>
00035 #include <deque>
00036 #include <cmath>
00037 
00038 #include "viennacl/forwards.h"
00039 
00040 namespace viennacl
00041 {
00042 
00043   namespace detail
00044   {
00045 
00046     // Calculate the bandwidth of a reordered matrix
00047     template <typename IndexT, typename ValueT>
00048     IndexT calc_reordered_bw(std::vector< std::map<IndexT, ValueT> > const & matrix,
00049                              std::vector<bool> & dof_assigned_to_node,
00050                              std::vector<IndexT> const & permutation)
00051     {
00052       IndexT bw = 0;
00053 
00054       for (vcl_size_t i = 0; i < permutation.size(); i++)
00055       {
00056         if (!dof_assigned_to_node[i])
00057           continue;
00058 
00059         IndexT min_index = static_cast<IndexT>(matrix.size());
00060         IndexT max_index = 0;
00061         for (typename std::map<IndexT, ValueT>::const_iterator it = matrix[i].begin(); it != matrix[i].end(); it++)
00062         {
00063           if (!dof_assigned_to_node[it->first])
00064             continue;
00065 
00066           if (permutation[it->first] > max_index)
00067             max_index = permutation[it->first];
00068           if (permutation[it->first] < min_index)
00069             min_index = permutation[it->first];
00070         }
00071         if (max_index > min_index)
00072           bw = std::max(bw, max_index - min_index);
00073       }
00074 
00075       return bw;
00076     }
00077 
00078 
00079 
00080     // function to calculate the increment of combination comb.
00081     // parameters:
00082     // comb: pointer to vector<int> of size m, m <= n
00083     //       1 <= comb[i] <= n for 0 <= i < m
00084     //       comb[i] < comb[i+1] for 0 <= i < m - 1
00085     //       comb represents an ordered selection of m values out of n
00086     // n: int
00087     //    total number of values out of which comb is taken as selection
00088     template <typename IndexT>
00089     bool comb_inc(std::vector<IndexT> & comb, vcl_size_t n)
00090     {
00091       IndexT m;
00092       IndexT k;
00093 
00094       m = static_cast<IndexT>(comb.size());
00095       // calculate k as highest possible index such that (*comb)[k-1] can be incremented
00096       k = m;
00097       while ( (k > 0) && ( ((k == m) && (comb[k-1] == static_cast<IndexT>(n)-1)) ||
00098                            ((k <  m) && (comb[k-1] == comb[k] - 1) )) )
00099       {
00100         k--;
00101       }
00102 
00103       if (k == 0) // no further increment of comb possible -> return false
00104         return false;
00105 
00106       comb[k-1] += 1;
00107 
00108       // and all higher index positions of comb are calculated just as directly following integer values
00109       // Example (1, 4, 7) -> (1, 5, 6) -> (1, 5, 7) -> (1, 6, 7) -> done   for n=7
00110       for (IndexT i = k; i < m; i++)
00111         comb[i] = comb[k-1] + (i - k);
00112       return true;
00113     }
00114 
00115 
00120     // node s
00121     template <typename MatrixT, typename IndexT>
00122     void generate_layering(MatrixT const & matrix,
00123                            std::vector< std::vector<IndexT> > & layer_list)
00124     {
00125       std::vector<bool> node_visited_already(matrix.size(), false);
00126 
00127       //
00128       // Step 1: Set root nodes to visited
00129       //
00130       for (vcl_size_t i=0; i<layer_list.size(); ++i)
00131       {
00132         for (typename std::vector<IndexT>::iterator it  = layer_list[i].begin();
00133                                                     it != layer_list[i].end();
00134                                                     it++)
00135           node_visited_already[*it] = true;
00136       }
00137 
00138       //
00139       // Step 2: Fill next layers
00140       //
00141       while (layer_list.back().size() > 0)
00142       {
00143         vcl_size_t layer_index = layer_list.size();  //parent nodes are at layer 0
00144         layer_list.push_back(std::vector<IndexT>());
00145 
00146         for (typename std::vector<IndexT>::iterator it  = layer_list[layer_index].begin();
00147                                                     it != layer_list[layer_index].end();
00148                                                     it++)
00149         {
00150           for (typename MatrixT::value_type::const_iterator it2  = matrix[*it].begin();
00151                                                             it2 != matrix[*it].end();
00152                                                             it2++)
00153           {
00154             if (it2->first == *it) continue;
00155             if (node_visited_already[it2->first]) continue;
00156 
00157             layer_list.back().push_back(it2->first);
00158             node_visited_already[it2->first] = true;
00159           }
00160         }
00161       }
00162 
00163       // remove last (empty) nodelist:
00164       layer_list.resize(layer_list.size()-1);
00165     }
00166 
00167 
00168     // function to generate a node layering as a tree structure rooted at node s
00169     template <typename MatrixType>
00170     void generate_layering(MatrixType const & matrix,
00171                            std::vector< std::vector<int> > & l,
00172                            int s)
00173     {
00174       vcl_size_t n = matrix.size();
00175       //std::vector< std::vector<int> > l;
00176       std::vector<bool> inr(n, false);
00177       std::vector<int> nlist;
00178 
00179       nlist.push_back(s);
00180       inr[s] = true;
00181       l.push_back(nlist);
00182 
00183       for (;;)
00184       {
00185           nlist.clear();
00186           for (std::vector<int>::iterator it  = l.back().begin();
00187                                           it != l.back().end();
00188                                           it++)
00189           {
00190               for (typename MatrixType::value_type::const_iterator it2  = matrix[*it].begin();
00191                                                          it2 != matrix[*it].end();
00192                                                          it2++)
00193               {
00194                   if (it2->first == *it) continue;
00195                   if (inr[it2->first]) continue;
00196 
00197                   nlist.push_back(it2->first);
00198                   inr[it2->first] = true;
00199               }
00200           }
00201 
00202           if (nlist.size() == 0)
00203               break;
00204 
00205           l.push_back(nlist);
00206       }
00207 
00208     }
00209 
00214     template <typename MatrixT, typename IndexT>
00215     void nodes_of_strongly_connected_component(MatrixT const & matrix,
00216                                                std::vector<IndexT> & node_list)
00217     {
00218       std::vector<bool> node_visited_already(matrix.size(), false);
00219       std::deque<IndexT> node_queue;
00220 
00221       //
00222       // Step 1: Push root nodes to queue:
00223       //
00224       for (typename std::vector<IndexT>::iterator it  = node_list.begin();
00225                                                   it != node_list.end();
00226                                                   it++)
00227       {
00228         node_queue.push_back(*it);
00229       }
00230       node_list.resize(0);
00231 
00232       //
00233       // Step 2: Fill with remaining nodes of strongly connected compontent
00234       //
00235       while (!node_queue.empty())
00236       {
00237         IndexT node_id = node_queue.front();
00238         node_queue.pop_front();
00239 
00240         if (!node_visited_already[node_id])
00241         {
00242           node_list.push_back(node_id);
00243           node_visited_already[node_id] = true;
00244 
00245           for (typename MatrixT::value_type::const_iterator it  = matrix[node_id].begin();
00246                                                             it != matrix[node_id].end();
00247                                                             it++)
00248           {
00249             IndexT neighbor_node_id = it->first;
00250             if (neighbor_node_id == node_id) continue;
00251             if (node_visited_already[neighbor_node_id]) continue;
00252 
00253             node_queue.push_back(neighbor_node_id);
00254           }
00255         }
00256       }
00257 
00258     }
00259 
00260 
00261     // comparison function for comparing two vector<int> values by their
00262     // [1]-element
00263     inline bool cuthill_mckee_comp_func(std::vector<int> const & a,
00264                                         std::vector<int> const & b)
00265     {
00266       return (a[1] < b[1]);
00267     }
00268 
00269     template <typename IndexT>
00270     bool cuthill_mckee_comp_func_pair(std::pair<IndexT, IndexT> const & a,
00271                                       std::pair<IndexT, IndexT> const & b)
00272     {
00273         return (a.second < b.second);
00274     }
00275 
00286     template <typename IndexT, typename ValueT>
00287     vcl_size_t cuthill_mckee_on_strongly_connected_component(std::vector< std::map<IndexT, ValueT> > const & matrix,
00288                                                               std::deque<IndexT> & node_assignment_queue,
00289                                                               std::vector<bool>  & dof_assigned_to_node,
00290                                                               std::vector<IndexT> & permutation,
00291                                                               vcl_size_t current_dof)
00292     {
00293       typedef std::pair<IndexT, IndexT> NodeIdDegreePair; //first member is the node ID, second member is the node degree
00294 
00295       std::vector< NodeIdDegreePair > local_neighbor_nodes(matrix.size());
00296 
00297       while (!node_assignment_queue.empty())
00298       {
00299         // Grab first node from queue
00300         vcl_size_t node_id = node_assignment_queue.front();
00301         node_assignment_queue.pop_front();
00302 
00303         // Assign dof if a new dof hasn't been assigned yet
00304         if (!dof_assigned_to_node[node_id])
00305         {
00306           permutation[node_id] = static_cast<IndexT>(current_dof);  //TODO: Invert this!
00307           ++current_dof;
00308           dof_assigned_to_node[node_id] = true;
00309 
00310           //
00311           // Get all neighbors of that node:
00312           //
00313           vcl_size_t num_neighbors = 0;
00314           for (typename std::map<IndexT, ValueT>::const_iterator neighbor_it  = matrix[node_id].begin();
00315                                                                  neighbor_it != matrix[node_id].end();
00316                                                                ++neighbor_it)
00317           {
00318             if (!dof_assigned_to_node[neighbor_it->first])
00319             {
00320               local_neighbor_nodes[num_neighbors] = NodeIdDegreePair(neighbor_it->first, static_cast<IndexT>(matrix[neighbor_it->first].size()));
00321               ++num_neighbors;
00322             }
00323           }
00324 
00325           // Sort neighbors by increasing node degree
00326           std::sort(local_neighbor_nodes.begin(), local_neighbor_nodes.begin() + num_neighbors, detail::cuthill_mckee_comp_func_pair<IndexT>);
00327 
00328           // Push neighbors to queue
00329           for (vcl_size_t i=0; i<num_neighbors; ++i)
00330             node_assignment_queue.push_back(local_neighbor_nodes[i].first);
00331 
00332         } // if node doesn't have a new dof yet
00333 
00334       } // while nodes in queue
00335 
00336       return current_dof;
00337 
00338     }
00339 
00340   } //namespace detail
00341 
00342   //
00343   // Part 1: The original Cuthill-McKee algorithm
00344   //
00345 
00347   struct cuthill_mckee_tag {};
00348 
00363   template <typename IndexT, typename ValueT>
00364   std::vector<IndexT> reorder(std::vector< std::map<IndexT, ValueT> > const & matrix, cuthill_mckee_tag)
00365   {
00366     std::vector<IndexT> permutation(matrix.size());
00367     std::vector<bool>   dof_assigned_to_node(matrix.size(), false);   //flag vector indicating whether node i has received a new dof
00368     std::deque<IndexT>  node_assignment_queue;
00369 
00370     vcl_size_t current_dof = 0;  //the dof to be assigned
00371 
00372     while (current_dof < matrix.size()) //outer loop for each strongly connected component (there may be more than one)
00373     {
00374       //
00375       // preprocessing: Determine node degrees for nodes which have not been assigned
00376       //
00377       vcl_size_t current_min_degree = matrix.size();
00378       vcl_size_t node_with_minimum_degree = 0;
00379       bool found_unassigned_node = false;
00380       for (vcl_size_t i=0; i<matrix.size(); ++i)
00381       {
00382         if (!dof_assigned_to_node[i])
00383         {
00384           if (matrix[i].size() == 1)  //This is an isolated node, so assign DOF right away
00385           {
00386             permutation[i] = static_cast<IndexT>(current_dof);
00387             dof_assigned_to_node[i] = true;
00388             ++current_dof;
00389             continue;
00390           }
00391 
00392           if (!found_unassigned_node) //initialize minimum degree on first node without new dof
00393           {
00394             current_min_degree = matrix[i].size();
00395             node_with_minimum_degree = i;
00396             found_unassigned_node = true;
00397           }
00398 
00399           if (matrix[i].size() < current_min_degree) //found a node with smaller degree
00400           {
00401             current_min_degree = matrix[i].size();
00402             node_with_minimum_degree = i;
00403           }
00404         }
00405       }
00406 
00407       //
00408       // Stage 2: Distribute dofs on this closely connected (sub-)graph in a breath-first manner using one root node
00409       //
00410       if (found_unassigned_node) // there's work to be done
00411       {
00412         node_assignment_queue.push_back(static_cast<IndexT>(node_with_minimum_degree));
00413         current_dof = detail::cuthill_mckee_on_strongly_connected_component(matrix, node_assignment_queue, dof_assigned_to_node, permutation, current_dof);
00414       }
00415     }
00416 
00417     return permutation;
00418   }
00419 
00420 
00421   //
00422   // Part 2: Advanced Cuthill McKee
00423   //
00424 
00426   class advanced_cuthill_mckee_tag
00427   {
00428     public:
00447       advanced_cuthill_mckee_tag(double a = 0.0, vcl_size_t gmax = 1) : starting_node_param_(a), max_root_nodes_(gmax) {}
00448 
00449       double starting_node_param() const { return starting_node_param_;}
00450       void starting_node_param(double a) { if (a >= 0) starting_node_param_ = a; }
00451 
00452       vcl_size_t max_root_nodes() const { return max_root_nodes_; }
00453       void max_root_nodes(vcl_size_t gmax) { max_root_nodes_ = gmax; }
00454 
00455     private:
00456       double starting_node_param_;
00457       vcl_size_t max_root_nodes_;
00458   };
00459 
00460 
00461 
00470   template <typename IndexT, typename ValueT>
00471   std::vector<IndexT> reorder(std::vector< std::map<IndexT, ValueT> > const & matrix,
00472                               advanced_cuthill_mckee_tag const & tag)
00473   {
00474     vcl_size_t n = matrix.size();
00475     double a = tag.starting_node_param();
00476     vcl_size_t gmax = tag.max_root_nodes();
00477     std::vector<IndexT> permutation(n);
00478     std::vector<bool>   dof_assigned_to_node(n, false);
00479     std::vector<IndexT> nodes_in_strongly_connected_component;
00480     std::vector<IndexT> parent_nodes;
00481     vcl_size_t deg_min;
00482     vcl_size_t deg_max;
00483     vcl_size_t deg_a;
00484     vcl_size_t deg;
00485     std::vector<IndexT> comb;
00486 
00487     nodes_in_strongly_connected_component.reserve(n);
00488     parent_nodes.reserve(n);
00489     comb.reserve(n);
00490 
00491     vcl_size_t current_dof = 0;
00492 
00493     while (current_dof < matrix.size()) // for all strongly connected components
00494     {
00495       // get all nodes of the strongly connected component:
00496       nodes_in_strongly_connected_component.resize(0);
00497       for (vcl_size_t i = 0; i < n; i++)
00498       {
00499         if (!dof_assigned_to_node[i])
00500         {
00501           nodes_in_strongly_connected_component.push_back(static_cast<IndexT>(i));
00502           detail::nodes_of_strongly_connected_component(matrix, nodes_in_strongly_connected_component);
00503           break;
00504         }
00505       }
00506 
00507       // determine minimum and maximum node degree
00508       deg_min = 0;
00509       deg_max = 0;
00510       for (typename std::vector<IndexT>::iterator it  = nodes_in_strongly_connected_component.begin();
00511                                                   it != nodes_in_strongly_connected_component.end();
00512                                                   it++)
00513       {
00514         deg = matrix[*it].size();
00515         if (deg_min == 0 || deg < deg_min)
00516           deg_min = deg;
00517         if (deg_max == 0 || deg > deg_max)
00518           deg_max = deg;
00519       }
00520       deg_a = deg_min + static_cast<vcl_size_t>(a * (deg_max - deg_min));
00521 
00522       // fill array of parent nodes:
00523       parent_nodes.resize(0);
00524       for (typename std::vector<IndexT>::iterator it  = nodes_in_strongly_connected_component.begin();
00525                                                   it != nodes_in_strongly_connected_component.end();
00526                                                   it++)
00527       {
00528         if (matrix[*it].size() <= deg_a)
00529           parent_nodes.push_back(*it);
00530       }
00531 
00532       //
00533       // backup current state in order to restore for every new combination of parent nodes below
00534       //
00535       std::vector<bool> dof_assigned_to_node_backup = dof_assigned_to_node;
00536       std::vector<bool> dof_assigned_to_node_best;
00537 
00538       std::vector<IndexT> permutation_backup = permutation;
00539       std::vector<IndexT> permutation_best = permutation;
00540 
00541       vcl_size_t current_dof_backup = current_dof;
00542 
00543       vcl_size_t g = 1;
00544       comb.resize(1);
00545       comb[0] = 0;
00546 
00547       IndexT bw_best = 0;
00548 
00549       //
00550       // Loop over all combinations of g <= gmax root nodes
00551       //
00552 
00553       for (;;)
00554       {
00555         dof_assigned_to_node = dof_assigned_to_node_backup;
00556         permutation          = permutation_backup;
00557         current_dof          = current_dof_backup;
00558 
00559         std::deque<IndexT>  node_queue;
00560 
00561         // add the selected root nodes according to actual combination comb to q
00562         for (typename std::vector<IndexT>::iterator it = comb.begin(); it != comb.end(); it++)
00563           node_queue.push_back(parent_nodes[*it]);
00564 
00565         current_dof = detail::cuthill_mckee_on_strongly_connected_component(matrix, node_queue, dof_assigned_to_node, permutation, current_dof);
00566 
00567         // calculate resulting bandwith for root node combination
00568         // comb for current numbered component of the node graph
00569         IndexT bw = detail::calc_reordered_bw(matrix, dof_assigned_to_node, permutation);
00570 
00571         // remember best ordering:
00572         if (bw_best == 0 || bw < bw_best)
00573         {
00574           permutation_best = permutation;
00575           bw_best = bw;
00576           dof_assigned_to_node_best = dof_assigned_to_node;
00577         }
00578 
00579         // calculate next combination comb, if not existing
00580         // increment g if g stays <= gmax, or else terminate loop
00581         if (!detail::comb_inc(comb, parent_nodes.size()))
00582         {
00583           ++g;
00584           if ( (gmax > 0 && g > gmax) || g > parent_nodes.size())
00585             break;
00586 
00587           comb.resize(g);
00588           for (vcl_size_t i = 0; i < g; i++)
00589             comb[i] = static_cast<IndexT>(i);
00590         }
00591       }
00592 
00593       //
00594       // restore best permutation
00595       //
00596       permutation = permutation_best;
00597       dof_assigned_to_node = dof_assigned_to_node_best;
00598 
00599     }
00600 
00601     return permutation;
00602   }
00603 
00604 
00605 } //namespace viennacl
00606 
00607 
00608 #endif