MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_NotayAggregationFactory_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_NOTAYAGGREGATIONFACTORY_DECL_HPP_
11 #define MUELU_NOTAYAGGREGATIONFACTORY_DECL_HPP_
12 
13 #include "MueLu_ConfigDefs.hpp"
14 
15 #include <Xpetra_Map_fwd.hpp>
16 #include <Xpetra_Vector_fwd.hpp>
17 
18 #include <Xpetra_Matrix_fwd.hpp>
19 
20 #include "MueLu_LWGraph_fwd.hpp"
22 #include "MueLu_Exceptions.hpp"
24 
26 
27 #include "MueLu_Level_fwd.hpp"
28 #include "MueLu_Aggregates_fwd.hpp"
29 #include "MueLu_Utilities_fwd.hpp"
30 
31 namespace MueLu {
32 
33 template <class Scalar = DefaultScalar,
36  class Node = DefaultNode>
38 #undef MUELU_NOTAYAGGREGATIONFACTORY_SHORT
39 #include "MueLu_UseShortNames.hpp"
40 
41  public:
43 
44  using local_matrix_type = typename Matrix::local_matrix_type;
45  using device_type = typename local_matrix_type::device_type;
46  using execution_space = typename device_type::execution_space;
48  using impl_scalar_type = typename Kokkos::ArithTraits<Scalar>::val_type;
49  using row_sum_type = typename Kokkos::View<impl_scalar_type*, Kokkos::LayoutLeft, device_type>;
51 
53 
54 
57 
60 
62 
64 
66 
67 
68  // Options shared by all aggregation algorithms
69 
71 
72 
73  void DeclareInput(Level& currentLevel) const;
74 
76 
78 
79 
81  void Build(Level& currentLevel) const;
82 
85  const RCP<const Matrix>& A,
86  const ArrayView<const LO>& orderingVector,
87  const magnitude_type kappa,
88  Aggregates& aggregates,
89  std::vector<unsigned>& aggStat,
90  LO& numNonAggregatedNodes,
91  LO& numDirichletNodes) const;
92 
95  const RCP<const Matrix>& A,
96  const Teuchos::ArrayView<const LO>& orderingVector,
97  const local_matrix_type& coarseA,
98  const magnitude_type kappa,
99  const row_sum_type& rowSum,
100  std::vector<LO>& localAggStat,
101  Array<LO>& localVertex2AggID,
102  LO& numLocalAggregates,
103  LO& numNonAggregatedNodes) const;
104 
105  void BuildOnRankLocalMatrix(const local_matrix_type& localA,
106  local_matrix_type& onRankA) const;
107 
109  void BuildIntermediateProlongator(const LO numRows,
110  const LO numDirichletNodes,
111  const LO numLocalAggregates,
112  const ArrayView<const LO>& localVertex2AggID,
113  local_matrix_type& intermediateP) const;
114 
116  void BuildCoarseLocalMatrix(const local_matrix_type& intermediateP,
117  local_matrix_type& coarseA) const;
118 
120  void localSpGEMM(const local_matrix_type& A,
121  const local_matrix_type& B,
122  const std::string matrixLabel,
123  local_matrix_type& C) const;
124 
126 
127  private:
128 }; // class NotayAggregationFactory
129 
130 } // namespace MueLu
131 
132 #define MUELU_NOTAYAGGREGATIONFACTORY_SHORT
133 #endif /* MUELU_NOTAYAGGREGATIONFACTORY_DECL_HPP_ */
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
MueLu::DefaultLocalOrdinal LocalOrdinal
void Build(Level &currentLevel) const
Build aggregates.
typename Kokkos::ArithTraits< Scalar >::val_type impl_scalar_type
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Container class for aggregation information.
typename Matrix::local_matrix_type local_matrix_type
typename Teuchos::ScalarTraits< Scalar >::magnitudeType magnitude_type
void BuildInitialAggregates(const Teuchos::ParameterList &params, const RCP< const Matrix > &A, const ArrayView< const LO > &orderingVector, const magnitude_type kappa, Aggregates &aggregates, std::vector< unsigned > &aggStat, LO &numNonAggregatedNodes, LO &numDirichletNodes) const
Initial aggregation phase.
typename local_matrix_type::device_type device_type
LocalOrdinal LO
void DeclareInput(Level &currentLevel) const
Input.
void BuildIntermediateProlongator(const LO numRows, const LO numDirichletNodes, const LO numLocalAggregates, const ArrayView< const LO > &localVertex2AggID, local_matrix_type &intermediateP) const
Construction of a local prolongator with values equal to 1.0.
MueLu::DefaultNode Node
void BuildFurtherAggregates(const Teuchos::ParameterList &params, const RCP< const Matrix > &A, const Teuchos::ArrayView< const LO > &orderingVector, const local_matrix_type &coarseA, const magnitude_type kappa, const row_sum_type &rowSum, std::vector< LO > &localAggStat, Array< LO > &localVertex2AggID, LO &numLocalAggregates, LO &numNonAggregatedNodes) const
Further aggregation phase increases coarsening rate by a factor of ~2 per iteration.
void BuildCoarseLocalMatrix(const local_matrix_type &intermediateP, local_matrix_type &coarseA) const
Implementation of a local Galerkin projection called inside BuildFurtherAggregates.
MueLu::DefaultScalar Scalar
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
typename device_type::execution_space execution_space
void BuildOnRankLocalMatrix(const local_matrix_type &localA, local_matrix_type &onRankA) const
void localSpGEMM(const local_matrix_type &A, const local_matrix_type &B, const std::string matrixLabel, local_matrix_type &C) const
Wrapper for kokkos-kernels&#39; spgemm that takes in CrsMatrix.
Base class for factories that use one level (currentLevel).
typename Kokkos::View< impl_scalar_type *, Kokkos::LayoutLeft, device_type > row_sum_type