MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_InterfaceAggregationAlgorithm_def.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/*
11 * MueLu_OnePtAggregationAlgorithm_def.hpp
12 *
13 * Created on: Sep 18, 2012
14 * Author: Tobias Wiesner
15 */
16
17#ifndef MUELU_INTERFACEAGGREGATIONALGORITHM_DEF_HPP_
18#define MUELU_INTERFACEAGGREGATIONALGORITHM_DEF_HPP_
19
20#include <Teuchos_Comm.hpp>
21#include <Teuchos_CommHelpers.hpp>
22
23#include <Xpetra_Vector.hpp>
24
26
27#include "MueLu_LWGraph.hpp"
28#include "MueLu_Aggregates.hpp"
29#include "MueLu_Exceptions.hpp"
30#include "MueLu_Monitor.hpp"
31
32namespace MueLu {
33
34template <class LocalOrdinal, class GlobalOrdinal, class Node>
37
38template <class LocalOrdinal, class GlobalOrdinal, class Node>
39void InterfaceAggregationAlgorithm<LocalOrdinal, GlobalOrdinal, Node>::BuildAggregatesNonKokkos(Teuchos::ParameterList const& /* params */, LWGraph const& graph, Aggregates& aggregates, typename AggregationAlgorithmBase<LocalOrdinal, GlobalOrdinal, Node>::AggStatHostType& aggStat, LO& numNonAggregatedNodes) const {
40 Monitor m(*this, "BuildAggregatesNonKokkos");
41
42 const LocalOrdinal nRows = graph.GetNodeNumVertices();
43 const int myRank = graph.GetComm()->getRank();
44
45 // vertex ids for output
46 Teuchos::ArrayRCP<LocalOrdinal> vertex2AggId = aggregates.GetVertex2AggId()->getDataNonConst(0);
47 Teuchos::ArrayRCP<LocalOrdinal> procWinner = aggregates.GetProcWinner()->getDataNonConst(0);
48
49 // some internal variables
50 LocalOrdinal numLocalAggregates = aggregates.GetNumAggregates(); // number of local aggregates on current proc
51
52 // main loop over all local rows of graph(A)
53 for (int iNode1 = 0; iNode1 < nRows; ++iNode1) {
54 if (aggStat[iNode1] == INTERFACE) {
55 aggregates.SetIsRoot(iNode1); // mark iNode1 as root node for new aggregate 'agg'
56 int aggIndex = numLocalAggregates;
57 std::vector<int> aggList;
58 aggList.push_back(iNode1);
59 auto neighOfINode = graph.getNeighborVertices(iNode1);
60
61 for (int j = 0; j < neighOfINode.length; ++j) {
62 LO neigh = neighOfINode(j);
63 if (neigh != iNode1 && graph.isLocalNeighborVertex(neigh)) {
64 if (aggStat[neigh] != AGGREGATED && aggStat[neigh] != INTERFACE &&
65 aggStat[neigh] != IGNORED) {
66 aggList.push_back(neigh);
67 }
68 }
69 }
70
71 for (size_t k = 0; k < aggList.size(); k++) {
72 aggStat[aggList[k]] = AGGREGATED;
73 vertex2AggId[aggList[k]] = aggIndex;
74 procWinner[aggList[k]] = myRank;
75 }
76 ++numLocalAggregates;
77 numNonAggregatedNodes -= aggList.size();
78 }
79
80 } // end for
81
82 // update aggregate object
83 aggregates.SetNumAggregates(numLocalAggregates);
84}
85
86template <class LocalOrdinal, class GlobalOrdinal, class Node>
87void InterfaceAggregationAlgorithm<LocalOrdinal, GlobalOrdinal, Node>::BuildAggregates(Teuchos::ParameterList const& /* params */, LWGraph_kokkos const& graph, Aggregates& aggregates, typename AggregationAlgorithmBase<LocalOrdinal, GlobalOrdinal, Node>::AggStatType& aggStat, LO& numNonAggregatedNodes) const {
88 TEUCHOS_ASSERT(false);
89}
90
91} // namespace MueLu
92
93#endif /* MUELU_INTERFACEAGGREGATIONALGORITHM_DEF_HPP_ */
MueLu::DefaultLocalOrdinal LocalOrdinal
Container class for aggregation information.
KOKKOS_INLINE_FUNCTION LO GetNumAggregates() const
void SetIsRoot(LO i, bool value=true)
Set root node information.
const RCP< LOMultiVector > & GetVertex2AggId() const
Returns constant vector that maps local node IDs to local aggregates IDs.
const RCP< LOVector > & GetProcWinner() const
Returns constant vector that maps local node IDs to owning processor IDs.
void SetNumAggregates(LO nAggregates)
Set number of local aggregates on current processor.
Kokkos::View< unsigned *, typename LWGraphHostType::device_type > AggStatHostType
Kokkos::View< unsigned *, typename LWGraphType::device_type > AggStatType
void BuildAggregatesNonKokkos(Teuchos::ParameterList const &params, LWGraph const &graph, Aggregates &aggregates, typename AggregationAlgorithmBase< LocalOrdinal, GlobalOrdinal, Node >::AggStatHostType &aggStat, LO &numNonAggregatedNodes) const
Local aggregation.
void BuildAggregates(const Teuchos::ParameterList &params, const LWGraph_kokkos &graph, Aggregates &aggregates, typename AggregationAlgorithmBase< LocalOrdinal, GlobalOrdinal, Node >::AggStatType &aggStat, LO &numNonAggregatedNodes) const
InterfaceAggregationAlgorithm(RCP< const FactoryBase > const &graphFact=Teuchos::null)
Constructor.
KOKKOS_INLINE_FUNCTION size_type GetNodeNumVertices() const
Return number of graph vertices.
KOKKOS_INLINE_FUNCTION bool isLocalNeighborVertex(LO i) const
Return true if vertex with local id 'v' is on current process.
const RCP< const Teuchos::Comm< int > > GetComm() const
KOKKOS_INLINE_FUNCTION neighbor_vertices_type getNeighborVertices(LO i) const
Return the list of vertices adjacent to the vertex 'v'.
Lightweight MueLu representation of a compressed row storage graph.
Lightweight MueLu representation of a compressed row storage graph.
Timer to be used in non-factories.
Namespace for MueLu classes and methods.