10#ifndef MUELU_NOTAYAGGREGATIONFACTORY_DECL_HPP_
11#define MUELU_NOTAYAGGREGATIONFACTORY_DECL_HPP_
15#include <Xpetra_Map_fwd.hpp>
16#include <Xpetra_Vector_fwd.hpp>
18#include <Xpetra_Matrix_fwd.hpp>
38#undef MUELU_NOTAYAGGREGATIONFACTORY_SHORT
47 using magnitude_type =
typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
48#if KOKKOS_VERSION >= 40799
49 using impl_scalar_type =
typename KokkosKernels::ArithTraits<Scalar>::val_type;
53 using row_sum_type =
typename Kokkos::View<impl_scalar_type*, Kokkos::LayoutLeft, device_type>;
89 const RCP<const Matrix>& A,
90 const ArrayView<const LO>& orderingVector,
93 std::vector<unsigned>& aggStat,
94 LO& numNonAggregatedNodes,
95 LO& numDirichletNodes)
const;
99 const RCP<const Matrix>& A,
100 const Teuchos::ArrayView<const LO>& orderingVector,
104 std::vector<LO>& localAggStat,
105 Array<LO>& localVertex2AggID,
106 LO& numLocalAggregates,
107 LO& numNonAggregatedNodes)
const;
114 const LO numDirichletNodes,
115 const LO numLocalAggregates,
116 const ArrayView<const LO>& localVertex2AggID,
126 const std::string matrixLabel,
136#define MUELU_NOTAYAGGREGATIONFACTORY_SHORT
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Container class for aggregation information.
Class that holds all level-specific information.
typename Kokkos::View< impl_scalar_type *, Kokkos::LayoutLeft, device_type > row_sum_type
void BuildInitialAggregates(const Teuchos::ParameterList ¶ms, 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.
virtual ~NotayAggregationFactory()
Destructor.
void BuildFurtherAggregates(const Teuchos::ParameterList ¶ms, 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.
typename Teuchos::ScalarTraits< Scalar >::magnitudeType magnitude_type
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.
void BuildCoarseLocalMatrix(const local_matrix_type &intermediateP, local_matrix_type &coarseA) const
Implementation of a local Galerkin projection called inside BuildFurtherAggregates.
typename Matrix::local_matrix_device_type local_matrix_type
void BuildOnRankLocalMatrix(const local_matrix_type &localA, local_matrix_type &onRankA) const
typename device_type::execution_space execution_space
typename local_matrix_type::device_type device_type
typename Kokkos::ArithTraits< Scalar >::val_type impl_scalar_type
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
NotayAggregationFactory()
Constructor.
void DeclareInput(Level ¤tLevel) const
Input.
void Build(Level ¤tLevel) const
Build aggregates.
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' spgemm that takes in CrsMatrix.
Base class for factories that use one level (currentLevel).
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar