10 #ifndef MUELU_AGGREGATIONSTRUCTUREDALGORITHM_DECL_HPP_
11 #define MUELU_AGGREGATIONSTRUCTUREDALGORITHM_DECL_HPP_
23 #include "MueLu_LWGraph.hpp"
24 #include "MueLu_LWGraph_kokkos.hpp"
49 #undef MUELU_AGGREGATIONSTRUCTUREDALGORITHM_SHORT
55 using size_type =
typename local_graph_type::size_type;
59 using memory_space =
typename local_graph_type::device_type::memory_space;
61 using LOVectorView = decltype(std::declval<LOVector>().getDeviceLocalView(Xpetra::Access::ReadWrite));
83 LO& numNonAggregatedNodes)
const;
97 LO& numNonAggregatedNodes)
const;
103 const LO dofsPerNode,
107 std::string
description()
const {
return "Aggretation: structured algorithm"; }
118 Kokkos::View<unsigned*, device_type> aggStat,
122 KOKKOS_INLINE_FUNCTION
138 const LO numGhostedNodes,
const LO dofsPerNode,
143 KOKKOS_INLINE_FUNCTION
158 const LO dofsPerNode,
159 const int numInterpolationPoints,
const LO numLocalRows,
163 KOKKOS_INLINE_FUNCTION
164 void operator()(
const LO rowIdx,
GO& update,
const bool final)
const;
181 const int numDimensions,
182 const LO numGhostedNodes,
const LO dofsPerNode,
183 const int numInterpolationPoints,
189 KOKKOS_INLINE_FUNCTION
196 const LO dofsPerNode,
const int numInterpolationPoints,
201 const LO dofsPerNode,
const int numInterpolationPoints,
208 #define MUELU_AGGREGATIONSTRUCTUREDALGORITHM_SHORT
non_const_row_map_type rowPtr_
Kokkos::View< unsigned *, typename LWGraphHostType::device_type > AggStatHostType
constLOTupleView ghostedNodesPerDir_
MueLu::DefaultLocalOrdinal LocalOrdinal
Lightweight MueLu representation of a compressed row storage graph.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
void ComputeGraphDataConstant(const LWGraph &graph, RCP< IndexManager > &geoData, const LO dofsPerNode, const int numInterpolationPoints, ArrayRCP< size_t > &nnzOnRow, Array< size_t > &rowPtr, Array< LO > &colIndex) const
Container class for aggregation information.
typename std::conditional< OnHost, typename local_graph_device_type::HostMirror, local_graph_device_type >::type local_graph_type
fillAggregatesFunctor(RCP< IndexManager_kokkos > geoData, const int myRank, Kokkos::View< unsigned *, device_type > aggStat, LOVectorView vertex2AggID, LOVectorView procWinner)
KOKKOS_INLINE_FUNCTION void operator()(const LO nodeIdx) const
constIntTupleView coarseRate_
decltype(std::declval< LOVector >().getDeviceLocalView(Xpetra::Access::ReadWrite)) LOVectorView
Algorithm for coarsening a graph with structured aggregation.
Pure virtual base class for all MueLu aggregation algorithms.
std::string description() const
Return a simple one-line description of this object.
KOKKOS_INLINE_FUNCTION void operator()(const LO nodeIdx) const
typename Kokkos::View< const LO[3], device_type > constLOTupleView
constIntTupleView endRate_
const int numGhostedNodes_
computeGraphRowPtrFunctor(RCP< IndexManager_kokkos > geoData, const LO dofsPerNode, const int numInterpolationPoints, const LO numLocalRows, constIntTupleView coarseRate, constLOTupleView lFineNodesPerDir, non_const_row_map_type rowPtr)
const int numInterpolationPoints_
typename local_graph_type::device_type::memory_space memory_space
typename local_graph_type::device_type::execution_space execution_space
constIntTupleView coarseRate_
constLOTupleView lFineNodesPerDir_
computeGraphDataLinearFunctor(RCP< IndexManager_kokkos > geoData, const int numDimensions, const LO numGhostedNodes, const LO dofsPerNode, const int numInterpolationPoints, constIntTupleView coarseRate, constIntTupleView endRate, constLOTupleView lFineNodesPerDir, constLOTupleView ghostedNodesPerDir, non_const_row_map_type rowPtr, entries_type colIndex)
virtual ~AggregationStructuredAlgorithm()
Destructor.
IndexManager_kokkos geoData_
const int numInterpolationPoints_
MueLu::DefaultGlobalOrdinal GlobalOrdinal
non_const_row_map_type rowPtr_
void BuildGraphOnHost(const LWGraph &graph, RCP< IndexManager > &geoData, const LO dofsPerNode, RCP< CrsGraph > &myGraph, RCP< const Map > &coarseCoordinatesFineMap, RCP< const Map > &coarseCoordinatesMap) const
Local aggregation.
IndexManager_kokkos geoData_
AggregationStructuredAlgorithm(const RCP< const FactoryBase > &=Teuchos::null)
Constructor.
void BuildAggregates(const Teuchos::ParameterList ¶ms, const LWGraph_kokkos &graph, Aggregates &aggregates, typename AggregationAlgorithmBase< LocalOrdinal, GlobalOrdinal, Node >::AggStatType &aggStat, LO &numNonAggregatedNodes) const
Build aggregates object.
typename local_graph_type::row_map_type::non_const_type non_const_row_map_type
IndexManager_kokkos geoData_
Kokkos::View< unsigned *, device_type > aggStat_
IndexManager_kokkos geoData_
KOKKOS_INLINE_FUNCTION void operator()(const LO rowIdx, GO &update, const bool final) const
typename LWGraph_kokkos::local_graph_type local_graph_type
KOKKOS_INLINE_FUNCTION void operator()(const LO nodeIdx, LO &lNumAggregatedNodes) const
constIntTupleView coarseRate_
void BuildAggregatesNonKokkos(const Teuchos::ParameterList ¶ms, const LWGraph &graph, Aggregates &aggregates, typename AggregationAlgorithmBase< LocalOrdinal, GlobalOrdinal, Node >::AggStatHostType &aggStat, LO &numNonAggregatedNodes) const
Local aggregation.
LOVectorView vertex2AggID_
const int numGhostedNodes_
typename local_graph_type::device_type device_type
typename local_graph_type::entries_type entries_type
computeGraphDataConstantFunctor(RCP< IndexManager_kokkos > geoData, const LO numGhostedNodes, const LO dofsPerNode, constIntTupleView coarseRate, constIntTupleView endRate, constLOTupleView lFineNodesPerDir, non_const_row_map_type rowPtr, entries_type colIndex)
constIntTupleView endRate_
constLOTupleView lFineNodesPerDir_
Lightweight MueLu representation of a compressed row storage graph.
Container class for mesh layout and indices calculation.
void ComputeGraphDataLinear(const LWGraph &graph, RCP< IndexManager > &geoData, const LO dofsPerNode, const int numInterpolationPoints, ArrayRCP< size_t > &nnzOnRow, Array< size_t > &rowPtr, Array< LO > &colIndex) const
non_const_row_map_type rowPtr_
void BuildGraph(const LWGraph_kokkos &graph, RCP< IndexManager_kokkos > &geoData, const LO dofsPerNode, RCP< CrsGraph > &myGraph) const
Build a CrsGraph instead of aggregates.
typename Kokkos::View< const int[3], device_type > constIntTupleView
Kokkos::View< unsigned *, typename LWGraphType::device_type > AggStatType
constLOTupleView lFineNodesPerDir_
typename local_graph_type::size_type size_type