• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

LCMPartition.h

Go to the documentation of this file.
00001 //*****************************************************************//
00002 //    Albany 2.0:  Copyright 2012 Sandia Corporation               //
00003 //    This Software is released under the BSD license detailed     //
00004 //    in the file "license.txt" in the top-level Albany directory  //
00005 //*****************************************************************//
00006 
00007 // Define only if Zoltan is enabled
00008 #if defined (ALBANY_LCM) && defined(ALBANY_ZOLTAN)
00009 
00010 #if !defined(LCM_Partition_h)
00011 #define LCM_Partition_h
00012 
00013 #include <iostream>
00014 #include <iterator>
00015 #include <map>
00016 #include <set>
00017 #include <vector>
00018 
00019 #include <boost/shared_ptr.hpp>
00020 
00021 #include <zoltan_cpp.h>
00022 
00023 #include <Teuchos_RCP.hpp>
00024 #include <Teuchos_ArrayRCP.hpp>
00025 #include <Teuchos_ParameterList.hpp>
00026 
00027 #include <Albany_AbstractDiscretization.hpp>
00028 #include <Albany_DiscretizationFactory.hpp>
00029 #include <Albany_STKDiscretization.hpp>
00030 #include <Albany_SolverFactory.hpp>
00031 #include <Albany_Utils.hpp>
00032 
00033 #include <Intrepid_MiniTensor_Geometry.h>
00034 
00035 using Intrepid::Index;
00036 using Intrepid::Vector;
00037 
00038 namespace LCM {
00039 
00043   typedef std::vector<int>
00044   IDList;
00045 
00050   typedef std::map<int, IDList>
00051   AdjacencyMap;
00052 
00056   typedef std::map<int, double>
00057   ScalarMap;
00058 
00062   typedef std::map<int, Vector<double> >
00063   PointMap;
00064 
00068   namespace PARTITION {
00069 
00070     enum Scheme {
00071       UNKNOWN,
00072       RANDOM,
00073       GEOMETRIC,
00074       HYPERGRAPH,
00075       KMEANS,
00076       SEQUENTIAL,
00077       KDTREE};
00078 
00079   }
00080 
00081   //
00083   //
00084   class ConnectivityArray;
00085   class DualGraph;
00086   class ZoltanHyperGraph;
00087   struct KDTreeNode;
00088 
00096   struct ClusterCenter {
00097 
00098     Vector<double>
00099     position;
00100 
00101     Vector<double>
00102     weighted_centroid;
00103 
00104     Index
00105     count;
00106 
00107   };
00108 
00116   struct KDTreeNode {
00117     std::string
00118     name;
00119 
00120     boost::shared_ptr<KDTreeNode>
00121     parent;
00122 
00123     // Children
00124     boost::shared_ptr<KDTreeNode>
00125     left;
00126 
00127     boost::shared_ptr<KDTreeNode>
00128     right;
00129 
00130     // Bounding box of cell
00131     Vector<double>
00132     lower_corner;
00133 
00134     Vector<double>
00135     upper_corner;
00136 
00137     // Weighted centroid and count
00138     Vector<double>
00139     weighted_centroid;
00140 
00141     Index
00142     count;
00143 
00144     std::set<Index>
00145     cell_points;
00146 
00147     std::set<Index>
00148     candidate_centers;
00149 
00150     Index
00151     closest_center_to_midcell;
00152 
00153   };
00154 
00162   template<typename Node>
00163   class KDTree {
00164   public:
00165 
00166     KDTree(std::vector<Vector<double> > const & points,
00167         Index const number_centers);
00168 
00169     boost::shared_ptr<Node> &
00170     get_root() {return root_;};
00171 
00172   private:
00173 
00174     boost::shared_ptr<Node>
00175     root_;
00176 
00177   };
00178 
00184   template<typename Node>
00185   boost::shared_ptr<Node>
00186   BuildKDTree(std::vector< Vector<double> > const & points);
00187 
00193   template<typename Node>
00194   boost::shared_ptr<Node>
00195   CreateKDTreeNode(
00196       std::string const & name,
00197       boost::shared_ptr<Node> parent,
00198       std::vector< Vector<double> > const & points,
00199       std::set<Index> const & points_indices);
00200 
00205   template<typename Node, typename Visitor>
00206   void
00207   VisitTreeNode(Node & node, Visitor const & visitor);
00208 
00212   template<typename Tree, typename Visitor>
00213   void
00214   TraverseTree(Tree & tree, Visitor const & visitor);
00215 
00219   template<typename Node>
00220   struct OutputVisitor {
00221     void
00222     operator()(Node const & node) const;
00223 
00224     bool
00225     pre_stop(Node const & node) const;
00226 
00227     bool
00228     post_stop(Node const & node) const;
00229   };
00230 
00234   template<typename Node, typename Center>
00235   struct FilterVisitor {
00236 
00237     std::vector<Vector<double> > & points;
00238 
00239     std::vector<Center> & centers;
00240 
00241     FilterVisitor(std::vector<Vector<double> > & p, std::vector<Center> & c);
00242 
00243     void
00244     operator()(Node const & node) const;
00245 
00246     bool
00247     pre_stop(Node const & node) const;
00248 
00249     bool
00250     post_stop(Node const & node) const;
00251   };
00252 
00257   class ConnectivityArray {
00258   public:
00259 
00263     ConnectivityArray();
00264 
00270     ConnectivityArray(
00271         std::string const & input_file,
00272         std::string const & output_file);
00273 
00277     Index
00278     GetNumberNodes() const;
00279 
00283     Index
00284     GetNumberElements() const;
00285 
00289     Index
00290     GetDimension() const;
00291 
00296     Intrepid::ELEMENT::Type
00297     GetType() const;
00298 
00303     Index
00304     GetNodesPerElement() const;
00305 
00309     PointMap
00310     GetNodeList() const;
00311 
00315     AdjacencyMap
00316     GetConnectivity() const;
00317 
00321     ScalarMap
00322     GetVolumes() const;
00323 
00327     double
00328     GetVolume() const;
00329 
00333     std::map<int, int>
00334     GetPartitions() const;
00335 
00339     ScalarMap
00340     GetPartitionVolumes() const;
00341 
00345     std::vector< Vector<double> >
00346     GetPartitionCentroids() const;
00347 
00351     PointMap
00352     GetCentroids() const;
00353 
00357     std::pair<Vector<double>, Vector<double> >
00358     BoundingBox() const;
00359 
00363     void
00364     SetTolerance(double tolerance);
00365 
00369     double
00370     GetTolerance() const;
00371 
00375     void
00376     SetCellSize(double requested_cell_size);
00377 
00381     double
00382     GetCellSize() const;
00383 
00387     void
00388     SetMaximumIterations(Index maximum_iterations);
00389 
00393     Index
00394     GetMaximumIterations() const;
00395 
00399     void
00400     SetInitializerScheme(PARTITION::Scheme initializer_scheme);
00401 
00405     PARTITION::Scheme
00406     GetInitializerScheme() const;
00407 
00411     void
00412     CheckNullVolume() const;
00413 
00419     std::vector< Vector<double> >
00420     CreateGrid();
00421 
00425     Vector<int>
00426     PointToIndex(Vector<double> const & point) const;
00427 
00431     bool
00432     IsInsideMesh(Vector<double> const & point) const;
00433 
00440     bool
00441     IsInsideMeshByElement(Vector<double> const & point) const;
00442 
00449     Index
00450     GetNumberPartitions(double const length_scale) const;
00451 
00455     Albany::AbstractDiscretization &
00456     GetDiscretization();
00457 
00463     std::map<int, int>
00464     PartitionByCenters(std::vector< Vector<double> > const & centers);
00465 
00473     std::map<int, int>
00474     Partition(
00475         const PARTITION::Scheme partition_scheme,
00476         double const length_scale);
00477 
00484     std::map<int, int>
00485     PartitionHyperGraph(double const length_scale);
00486 
00493     std::map<int, int>
00494     PartitionGeometric(double const length_scale);
00495 
00502     std::map<int, int>
00503     PartitionKMeans(double const length_scale);
00504 
00511     std::map<int, int>
00512     PartitionKDTree(double const length_scale);
00513 
00520     std::map<int, int>
00521     PartitionSequential(double const length_scale);
00522 
00530     std::map<int, int>
00531     PartitionRandom(double const length_scale);
00532 
00545     static int
00546     GetNumberGeometry(
00547         void* data,
00548         int* ierr);
00549 
00560     static int
00561     GetNumberOfObjects(
00562         void* data,
00563         int* ierr);
00564 
00602     static void
00603     GetObjectList(
00604         void* data,
00605         int sizeGID,
00606         int sizeLID,
00607         ZOLTAN_ID_PTR globalID,
00608         ZOLTAN_ID_PTR localID,
00609         int wgt_dim,
00610         float* obj_wgts,
00611         int* ierr);
00612 
00651     static void
00652     GetGeometry(
00653         void* data,
00654         int sizeGID,
00655         int sizeLID,
00656         int num_obj,
00657         ZOLTAN_ID_PTR globalID,
00658         ZOLTAN_ID_PTR localID,
00659         int num_dim,
00660         double* geom_vec,
00661         int* ierr);
00662 
00663   private:
00664 
00665     //
00666     // The type of elements in the mesh (assumed that all are of same type)
00667     //
00668     Intrepid::ELEMENT::Type
00669     type_;
00670 
00671     //
00672     // Node list
00673     //
00674     PointMap
00675     nodes_;
00676 
00677     //
00678     // Element - nodes connectivity
00679     //
00680     AdjacencyMap
00681     connectivity_;
00682 
00683     //
00684     // Space dimension
00685     //
00686     Index
00687     dimension_;
00688 
00689     //
00690     // Teuchos pointer to corresponding discretization
00691     //
00692     Teuchos::RCP<Albany::AbstractDiscretization>
00693     discretization_ptr_;
00694 
00695     //
00696     // Partitions if mesh is partitioned; otherwise empty
00697     //
00698     std::map<int, int>
00699     partitions_;
00700 
00701     //
00702     // Background grid of the domain for fast determination
00703     // of whether a point is inside the domain or not.
00704     //
00705     std::vector< std::vector< std::vector<bool> > >
00706     cells_;
00707 
00708     //
00709     // Size of background grid cell
00710     //
00711     Vector<double>
00712     cell_size_;
00713 
00714     //
00715     // Parameters for kmeans partitioning
00716     //
00717     double
00718     tolerance_;
00719 
00720     double
00721     requested_cell_size_;
00722 
00723     Index
00724     maximum_iterations_;
00725 
00726     //
00727     // Limits of the bounding box for coordinate array
00728     //
00729     Vector<double>
00730     lower_corner_;
00731 
00732     Vector<double>
00733     upper_corner_;
00734 
00735     //
00736     // Initializer scheme, if any.
00737     //
00738     PARTITION::Scheme
00739     initializer_scheme_;
00740 
00741   };
00742 
00746   class DualGraph {
00747 
00748   public:
00749 
00753     DualGraph();
00754 
00758     DualGraph(ConnectivityArray const & connectivity_array);
00759 
00763     int
00764     GetNumberVertices() const;
00765 
00769     int
00770     GetNumberEdges() const;
00771 
00776     void
00777     SetVertexWeights(ScalarMap & vertex_weights);
00778 
00782     ScalarMap
00783     GetVertexWeights() const;
00784 
00789     void
00790     SetGraph(AdjacencyMap & graph);
00791 
00795     AdjacencyMap
00796     GetGraph() const;
00797 
00801     AdjacencyMap
00802     GetEdgeList() const;
00803 
00807     int
00808     GetConnectedComponents(std::vector<int> & components) const;
00809 
00813     void
00814     Print() const;
00815 
00816   private:
00817 
00818     //
00819     // Given a connectivity array type, return local numbering of
00820     // proper faces
00821     //
00822     std::vector< std::vector<int> >
00823     GetFaceConnectivity(Intrepid::ELEMENT::Type const type) const;
00824 
00825   private:
00826 
00827     //
00828     // Number of edges in dual graph
00829     //
00830     int
00831     number_edges_;
00832 
00833     //
00834     // Graph data structure
00835     //
00836     AdjacencyMap
00837     graph_;
00838 
00839     //
00840     // Vertex weights
00841     //
00842     ScalarMap
00843     vertex_weights_;
00844 
00845   };
00846 
00858   class ZoltanHyperGraph {
00859 
00860   public:
00861 
00865     ZoltanHyperGraph();
00866 
00871     ZoltanHyperGraph(DualGraph const & dual_graph);
00872 
00876     int
00877     GetNumberVertices() const;
00878 
00883     void
00884     SetNumberHyperedges(int number_hyperedges);
00885 
00889     int
00890     GetNumberHyperedges() const;
00891 
00896     void
00897     SetGraph(AdjacencyMap & graph);
00898 
00902     AdjacencyMap
00903     GetGraph() const;
00904 
00909     void
00910     SetVertexWeights(ScalarMap & vertex_weights);
00911 
00915     ScalarMap
00916     GetVertexWeights() const;
00917 
00921     std::vector<ZOLTAN_ID_TYPE>
00922     GetEdgeIDs() const;
00923 
00928     std::vector<int>
00929     GetEdgePointers() const;
00930 
00934     std::vector<ZOLTAN_ID_TYPE>
00935     GetVertexIDs() const;
00936 
00947     static int
00948     GetNumberOfObjects(
00949         void* data,
00950         int* ierr);
00951 
00989     static void
00990     GetObjectList(
00991         void* data,
00992         int sizeGID,
00993         int sizeLID,
00994         ZOLTAN_ID_PTR globalID,
00995         ZOLTAN_ID_PTR localID,
00996         int wgt_dim,
00997         float* obj_wgts,
00998         int* ierr);
00999 
01024     static void
01025     GetHyperGraphSize(
01026         void* data,
01027         int* num_lists,
01028         int* num_pins,
01029         int* format,
01030         int* ierr);
01031 
01078     static void
01079     GetHyperGraph(
01080         void* data,
01081         int num_gid_entries,
01082         int num_vtx_edge,
01083         int num_pins,
01084         int format,
01085         ZOLTAN_ID_PTR
01086         vtxedge_GID,
01087         int* vtxedge_ptr,
01088         ZOLTAN_ID_PTR pin_GID,
01089         int* ierr);
01090 
01091   private:
01092 
01093     //
01094     // Number of vertices
01095     //
01096     int
01097     number_vertices_;
01098 
01099     //
01100     // Number of hyperedges
01101     //
01102     int
01103     number_hyperedges_;
01104 
01105     //
01106     // Graph data structure
01107     //
01108     AdjacencyMap
01109     graph_;
01110 
01111     //
01112     // Vertex weights
01113     //
01114     ScalarMap
01115     vertex_weights_;
01116 
01117   };
01118 
01124   std::istream &
01125   operator>>(
01126       std::istream & input_stream,
01127       ConnectivityArray & connectivity_array);
01128 
01134   std::ostream &
01135   operator<<(
01136       std::ostream & output_stream,
01137       ConnectivityArray const & connectivity_array);
01138 
01144   std::istream &
01145   operator>>(
01146       std::istream & input_stream,
01147       ZoltanHyperGraph & zoltan_hypergraph);
01148 
01154   std::ostream &
01155   operator<<(
01156       std::ostream & output_stream,
01157       ZoltanHyperGraph const & zoltan_hypergraph);
01158 
01159 } // namespace LCM
01160 
01161 #endif // #if !defined(LCM_Partition_h)
01162 
01163 #endif // #if defined (ALBANY_LCM) && defined(ALBANY_ZOLTAN)

Generated on Wed Mar 26 2014 18:36:39 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1