MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_AggregationStructuredAlgorithm_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_AGGREGATIONSTRUCTUREDALGORITHM_DECL_HPP_
11#define MUELU_AGGREGATIONSTRUCTUREDALGORITHM_DECL_HPP_
12
13#include "MueLu_ConfigDefs.hpp"
16
17#include "Xpetra_Vector.hpp"
18
23#include "MueLu_LWGraph.hpp"
24#include "MueLu_LWGraph_kokkos.hpp"
25
26namespace MueLu {
45template <class LocalOrdinal = DefaultLocalOrdinal,
47 class Node = DefaultNode>
48class AggregationStructuredAlgorithm : public MueLu::AggregationAlgorithmBase<LocalOrdinal, GlobalOrdinal, Node> {
49#undef MUELU_AGGREGATIONSTRUCTUREDALGORITHM_SHORT
51
52 public:
54 using non_const_row_map_type = typename local_graph_type::row_map_type::non_const_type;
55 using size_type = typename local_graph_type::size_type;
56 using entries_type = typename local_graph_type::entries_type;
57 using device_type = typename local_graph_type::device_type;
58 using execution_space = typename local_graph_type::device_type::execution_space;
59 using memory_space = typename local_graph_type::device_type::memory_space;
60
61 using LOVectorView = decltype(std::declval<LOVector>().getLocalViewDevice(Tpetra::Access::ReadWrite));
62 using constIntTupleView = typename Kokkos::View<const int[3], device_type>;
63 using constLOTupleView = typename Kokkos::View<const LO[3], device_type>;
65
66
68 AggregationStructuredAlgorithm(const RCP<const FactoryBase>& /* graphFact */ = Teuchos::null) {}
69
72
74
76
77
80 void BuildAggregatesNonKokkos(const Teuchos::ParameterList& params, const LWGraph& graph,
81 Aggregates& aggregates,
83 LO& numNonAggregatedNodes) const;
84
87 void BuildGraphOnHost(const LWGraph& graph, RCP<IndexManager>& geoData, const LO dofsPerNode,
88 RCP<CrsGraph>& myGraph, RCP<const Map>& coarseCoordinatesFineMap,
89 RCP<const Map>& coarseCoordinatesMap) const;
90
93 void BuildAggregates(const Teuchos::ParameterList& params,
94 const LWGraph_kokkos& graph,
95 Aggregates& aggregates,
97 LO& numNonAggregatedNodes) const;
98
101 void BuildGraph(const LWGraph_kokkos& graph,
102 RCP<IndexManager_kokkos>& geoData,
103 const LO dofsPerNode,
104 RCP<CrsGraph>& myGraph) const;
106
107 std::string description() const { return "Aggretation: structured algorithm"; }
108
111 const int myRank_;
112 Kokkos::View<unsigned*, device_type> aggStat_;
115
116 fillAggregatesFunctor(RCP<IndexManager_kokkos> geoData,
117 const int myRank,
118 Kokkos::View<unsigned*, device_type> aggStat,
119 LOVectorView vertex2AggID,
120 LOVectorView procWinner);
121
122 KOKKOS_INLINE_FUNCTION
123 void operator()(const LO nodeIdx, LO& lNumAggregatedNodes) const;
124
125 }; // struct fillAggregatesFunctor
126
130 const LO dofsPerNode_;
136
137 computeGraphDataConstantFunctor(RCP<IndexManager_kokkos> geoData,
138 const LO numGhostedNodes, const LO dofsPerNode,
139 constIntTupleView coarseRate, constIntTupleView endRate,
140 constLOTupleView lFineNodesPerDir,
141 non_const_row_map_type rowPtr, entries_type colIndex);
142
143 KOKKOS_INLINE_FUNCTION
144 void operator()(const LO nodeIdx) const;
145
146 }; // struct computeGraphDataConstantFunctor
147
150 const LO dofsPerNode_;
156
157 computeGraphRowPtrFunctor(RCP<IndexManager_kokkos> geoData,
158 const LO dofsPerNode,
159 const int numInterpolationPoints, const LO numLocalRows,
160 constIntTupleView coarseRate, constLOTupleView lFineNodesPerDir,
162
163 KOKKOS_INLINE_FUNCTION
164 void operator()(const LO rowIdx, GO& update, const bool final) const;
165 }; // struct computeGraphRowPtrFunctor
166
169 const int numDimensions_;
171 const LO dofsPerNode_;
179
180 computeGraphDataLinearFunctor(RCP<IndexManager_kokkos> geoData,
181 const int numDimensions,
182 const LO numGhostedNodes, const LO dofsPerNode,
183 const int numInterpolationPoints,
184 constIntTupleView coarseRate, constIntTupleView endRate,
185 constLOTupleView lFineNodesPerDir,
186 constLOTupleView ghostedNodesPerDir,
187 non_const_row_map_type rowPtr, entries_type colIndex);
188
189 KOKKOS_INLINE_FUNCTION
190 void operator()(const LO nodeIdx) const;
191
192 }; // struct computeGraphDataLinearFunctor
193
194 private:
195 void ComputeGraphDataConstant(const LWGraph& graph, RCP<IndexManager>& geoData,
196 const LO dofsPerNode, const int numInterpolationPoints,
197 ArrayRCP<size_t>& nnzOnRow, Array<size_t>& rowPtr,
198 Array<LO>& colIndex) const;
199
200 void ComputeGraphDataLinear(const LWGraph& graph, RCP<IndexManager>& geoData,
201 const LO dofsPerNode, const int numInterpolationPoints,
202 ArrayRCP<size_t>& nnzOnRow, Array<size_t>& rowPtr,
203 Array<LO>& colIndex) const;
204};
205
206} // namespace MueLu
207
208#define MUELU_AGGREGATIONSTRUCTUREDALGORITHM_SHORT
209#endif /* MUELU_AGGREGATIONSTRUCTUREDALGORITHM_DECL_HPP_ */
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Container class for aggregation information.
Pure virtual base class for all MueLu aggregation algorithms.
Kokkos::View< unsigned *, typename LWGraphHostType::device_type > AggStatHostType
Kokkos::View< unsigned *, typename LWGraphType::device_type > AggStatType
Algorithm for coarsening a graph with structured aggregation.
decltype(std::declval< LOVector >().getLocalViewDevice(Tpetra::Access::ReadWrite)) LOVectorView
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
void BuildGraph(const LWGraph_kokkos &graph, RCP< IndexManager_kokkos > &geoData, const LO dofsPerNode, RCP< CrsGraph > &myGraph) const
Build a CrsGraph instead of aggregates.
std::string description() const
Return a simple one-line description of this object.
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
typename Kokkos::View< const int[3], device_type > constIntTupleView
typename local_graph_type::device_type::execution_space execution_space
typename LWGraph_kokkos::local_graph_type local_graph_type
typename Kokkos::View< const LO[3], device_type > constLOTupleView
void BuildAggregatesNonKokkos(const Teuchos::ParameterList &params, const LWGraph &graph, Aggregates &aggregates, typename AggregationAlgorithmBase< LocalOrdinal, GlobalOrdinal, Node >::AggStatHostType &aggStat, LO &numNonAggregatedNodes) const
Local aggregation.
typename local_graph_type::row_map_type::non_const_type non_const_row_map_type
AggregationStructuredAlgorithm(const RCP< const FactoryBase > &=Teuchos::null)
Constructor.
void BuildAggregates(const Teuchos::ParameterList &params, const LWGraph_kokkos &graph, Aggregates &aggregates, typename AggregationAlgorithmBase< LocalOrdinal, GlobalOrdinal, Node >::AggStatType &aggStat, LO &numNonAggregatedNodes) const
Build aggregates object.
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.
typename local_graph_type::device_type::memory_space memory_space
Container class for mesh layout and indices calculation.
typename std::conditional< OnHost, typename local_graph_device_type::HostMirror, local_graph_device_type >::type local_graph_type
Lightweight MueLu representation of a compressed row storage graph.
Lightweight MueLu representation of a compressed row storage graph.
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
KOKKOS_INLINE_FUNCTION void operator()(const LO rowIdx, GO &update, const bool final) const
KOKKOS_INLINE_FUNCTION void operator()(const LO nodeIdx, LO &lNumAggregatedNodes) const