MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_BrickAggregationFactory_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // MueLu: A package for multigrid based preconditioning
4 //
5 // Copyright 2012 NTESS and the MueLu contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef MUELU_BRICKAGGREGATIONFACTORY_DECL_HPP_
11 #define MUELU_BRICKAGGREGATIONFACTORY_DECL_HPP_
12 
13 #include "MueLu_ConfigDefs.hpp"
14 
15 #include <Xpetra_Import_fwd.hpp>
17 #include <Xpetra_Map_fwd.hpp>
19 #include <Xpetra_Matrix_fwd.hpp>
22 
25 
26 #include "MueLu_LWGraph_fwd.hpp"
27 
28 #include "MueLu_LWGraph_fwd.hpp"
29 #include "MueLu_Level_fwd.hpp"
30 #include "MueLu_Aggregates_fwd.hpp"
31 #include "MueLu_Exceptions.hpp"
32 #include "MueLu_Utilities_fwd.hpp"
33 
42 namespace MueLu {
43 
44 template <class Scalar = DefaultScalar,
47  class Node = DefaultNode>
49 #undef MUELU_BRICKAGGREGATIONFACTORY_SHORT
50 #include "MueLu_UseShortNames.hpp"
51  private:
53 
54  // Comparator for doubles
55  // Generally, the coordinates for coarser levels would come out of averaging of fine level coordinates
56  // It is possible that the result of the averaging differs slightly between clusters, as we might have
57  // 3x2 and 2x2 cluster which would result in averaging 6 and 4 y-coordinates respectively, leading to
58  // slightly different results.
59  // Therefore, we hardcode a constant so that close points are considered the same.
60  class compare {
61  public:
62  bool operator()(const Scalar& x, const Scalar& y) const {
63  if (STS::magnitude(x - y) < 1e-14)
64  return false;
65  return STS::real(x) < STS::real(y);
66  }
67  };
68  typedef std::map<Scalar, GlobalOrdinal, compare> container;
69 
70  public:
72 
73 
76  : nDim_(-1)
77  , nx_(-1)
78  , ny_(-1)
79  , nz_(-1)
80  , bx_(-1)
81  , by_(-1)
82  , bz_(-1){};
83 
86 
88 
90 
91  // Options shared by all aggregation algorithms
92 
94 
95 
96  void DeclareInput(Level& currentLevel) const;
97 
99 
101 
102 
104  void Build(Level& currentLevel) const;
105 
107 
108  private:
109  void Setup(const RCP<const Teuchos::Comm<int> >& comm, const RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO> >& coords, const RCP<const Map>& map) const;
111 
112  void BuildGraph(Level& currentLevel, const RCP<Matrix>& A) const;
113 
114  bool isDirichlet(LocalOrdinal LID) const;
115  bool isRoot(LocalOrdinal LID) const;
118 
119  void getIJK(LocalOrdinal LID, int& i, int& j, int& k) const;
120  void getAggIJK(LocalOrdinal LID, int& i, int& j, int& k) const;
121 
122  mutable int nDim_;
125  mutable int nx_, ny_, nz_;
126  mutable int bx_, by_, bz_;
128  mutable int naggx_, naggy_, naggz_;
129 
130  mutable std::map<GlobalOrdinal, GlobalOrdinal> revMap_;
131 }; // class BrickAggregationFactory
132 
133 } // namespace MueLu
134 
135 #define MUELU_BRICKAGGREGATIONFACTORY_SHORT
136 #endif /* MUELU_BRICKAGGREGATIONFACTORY_DECL_HPP_ */
void getAggIJK(LocalOrdinal LID, int &i, int &j, int &k) const
ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::magnitudeType > x_
MueLu::DefaultLocalOrdinal LocalOrdinal
std::map< GlobalOrdinal, GlobalOrdinal > revMap_
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
GlobalOrdinal GO
static magnitudeType real(T a)
LocalOrdinal LO
void Build(Level &currentLevel) const
Build aggregates.
MueLu::DefaultNode Node
GlobalOrdinal getAggGID(LocalOrdinal LID) const
std::map< Scalar, GlobalOrdinal, compare > container
MueLu::DefaultScalar Scalar
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
void Setup(const RCP< const Teuchos::Comm< int > > &comm, const RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LO, GO, NO > > &coords, const RCP< const Map > &map) const
ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::magnitudeType > z_
static magnitudeType magnitude(T a)
void getIJK(LocalOrdinal LID, int &i, int &j, int &k) const
bool operator()(const Scalar &x, const Scalar &y) const
RCP< container > Construct1DMap(const RCP< const Teuchos::Comm< int > > &comm, const ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &x) const
Node NO
ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::magnitudeType > y_
void BuildGraph(Level &currentLevel, const RCP< Matrix > &A) const
Base class for factories that use one level (currentLevel).
GlobalOrdinal getRoot(LocalOrdinal LID) const
void DeclareInput(Level &currentLevel) const
Input.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.