MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_BlackBoxPFactory_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_BLACKBOXPFACTORY_DECL_HPP
11#define MUELU_BLACKBOXPFACTORY_DECL_HPP
12
13#include <Teuchos_SerialDenseVector.hpp>
14
15#include <Xpetra_MultiVector_fwd.hpp>
16#include <Xpetra_Matrix_fwd.hpp>
17
18#include "MueLu_ConfigDefs.hpp"
19#include "MueLu_PFactory.hpp"
21
22#include "MueLu_Level_fwd.hpp"
23
24namespace MueLuTests {
25// Forward declaration of friend tester class used to UnitTest BlackBoxPFactory
26template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
28} // namespace MueLuTests
29
30namespace MueLu {
31
80template <class Scalar = DefaultScalar,
83 class Node = DefaultNode>
84class BlackBoxPFactory : public PFactory {
85#undef MUELU_BLACKBOXPFACTORY_SHORT
87
88 public:
90
92
93
96
98 virtual ~BlackBoxPFactory() {}
100
101 RCP<const ParameterList> GetValidParameterList() const;
102
104
105
106 void DeclareInput(Level& fineLevel, Level& coarseLevel) const;
107
109
111
112
113 void Build(Level& fineLevel, Level& coarseLevel) const;
114 void BuildP(Level& fineLevel, Level& coarseLevel) const;
115
117
118 private:
119 struct NodesIDs {
120 // This small struct just carries basic data associated with coarse nodes that is needed
121 // to compute colMapP and to fillComplete P,
122
123 Array<GO> GIDs, coarseGIDs;
124 Array<int> PIDs;
125 Array<LO> LIDs;
126 std::vector<GO> colInds;
127 };
128
129 struct NodeID {
130 // This small struct is similar to the one above but only for one node.
131 // It is used to create a vector of NodeID that can easily be sorted
132
133 GO GID;
134 int PID;
136 };
137
138 void GetGeometricData(RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO> >& coordinates,
139 const Array<LO> coarseRate, const Array<GO> gFineNodesPerDir,
140 const Array<LO> lFineNodesPerDir, const LO BlkSize, Array<GO>& gIndices,
141 Array<LO>& myOffset, Array<bool>& ghostInterface, Array<LO>& endRate,
142 Array<GO>& gCoarseNodesPerDir, Array<LO>& lCoarseNodesPerDir,
143 Array<LO>& glCoarseNodesPerDir, Array<GO>& ghostGIDs,
144 Array<GO>& coarseNodesGIDs, Array<GO>& colGIDs, GO& gNumCoarseNodes,
145 LO& lNumCoarseNodes, ArrayRCP<Array<typename Teuchos::ScalarTraits<Scalar>::magnitudeType> > coarseNodes,
146 Array<int>& boundaryFlags, RCP<NodesIDs> ghostedCoarseNodes) const;
147
148 void ComputeLocalEntries(const RCP<const Matrix>& Aghost, const Array<LO> coarseRate,
149 const Array<LO> endRate, const LO BlkSize, const Array<LO> elemInds,
150 const Array<LO> lCoarseElementsPerDir,
151 const LO numDimensions, const Array<LO> lFineNodesPerDir,
152 const Array<GO> gFineNodesPerDir, const Array<GO> gIndices,
153 const Array<LO> lCoarseNodesPerDir, const Array<bool> ghostInterface,
154 const Array<int> elementFlags, const std::string stencilType,
155 const std::string blockStrategy, const Array<LO> elementNodesPerDir,
156 const LO numNodesInElement, const Array<GO> colGIDs,
157 Teuchos::SerialDenseMatrix<LO, SC>& Pi,
158 Teuchos::SerialDenseMatrix<LO, SC>& Pf,
159 Teuchos::SerialDenseMatrix<LO, SC>& Pe,
160 Array<LO>& dofType, Array<LO>& lDofInd) const;
161
162 void CollapseStencil(const int type, const int orientation, const int collapseFlags[3],
163 Array<SC>& stencil) const;
164
165 void FormatStencil(const LO BlkSize, const Array<bool> ghostInterface, const LO ie,
166 const LO je, const LO ke, const ArrayView<const SC> rowValues,
167 const Array<LO> elementNodesPerDir, const int collapseFlags[3],
168 const std::string stencilType, Array<SC>& stencil) const;
169
170 void GetNodeInfo(const LO ie, const LO je, const LO ke, const Array<LO> elementNodesPerDir,
171 int* type, LO& ind, int* orientation) const;
172
173 void sh_sort_permute(
174 const typename Teuchos::Array<LocalOrdinal>::iterator& first1,
175 const typename Teuchos::Array<LocalOrdinal>::iterator& last1,
176 const typename Teuchos::Array<LocalOrdinal>::iterator& first2,
177 const typename Teuchos::Array<LocalOrdinal>::iterator& last2) const;
178
179}; // class BlackBoxPFactory
180
181} // namespace MueLu
182
183#define MUELU_BLACKBOXPFACTORY_SHORT
184#endif // MUELU_BLACKBOXPFACTORY_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Prolongator factory performing geometric coarsening.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
void FormatStencil(const LO BlkSize, const Array< bool > ghostInterface, const LO ie, const LO je, const LO ke, const ArrayView< const SC > rowValues, const Array< LO > elementNodesPerDir, const int collapseFlags[3], const std::string stencilType, Array< SC > &stencil) const
void GetNodeInfo(const LO ie, const LO je, const LO ke, const Array< LO > elementNodesPerDir, int *type, LO &ind, int *orientation) const
void sh_sort_permute(const typename Teuchos::Array< LocalOrdinal >::iterator &first1, const typename Teuchos::Array< LocalOrdinal >::iterator &last1, const typename Teuchos::Array< LocalOrdinal >::iterator &first2, const typename Teuchos::Array< LocalOrdinal >::iterator &last2) const
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void ComputeLocalEntries(const RCP< const Matrix > &Aghost, const Array< LO > coarseRate, const Array< LO > endRate, const LO BlkSize, const Array< LO > elemInds, const Array< LO > lCoarseElementsPerDir, const LO numDimensions, const Array< LO > lFineNodesPerDir, const Array< GO > gFineNodesPerDir, const Array< GO > gIndices, const Array< LO > lCoarseNodesPerDir, const Array< bool > ghostInterface, const Array< int > elementFlags, const std::string stencilType, const std::string blockStrategy, const Array< LO > elementNodesPerDir, const LO numNodesInElement, const Array< GO > colGIDs, Teuchos::SerialDenseMatrix< LO, SC > &Pi, Teuchos::SerialDenseMatrix< LO, SC > &Pf, Teuchos::SerialDenseMatrix< LO, SC > &Pe, Array< LO > &dofType, Array< LO > &lDofInd) const
void CollapseStencil(const int type, const int orientation, const int collapseFlags[3], Array< SC > &stencil) const
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void GetGeometricData(RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LO, GO, NO > > &coordinates, const Array< LO > coarseRate, const Array< GO > gFineNodesPerDir, const Array< LO > lFineNodesPerDir, const LO BlkSize, Array< GO > &gIndices, Array< LO > &myOffset, Array< bool > &ghostInterface, Array< LO > &endRate, Array< GO > &gCoarseNodesPerDir, Array< LO > &lCoarseNodesPerDir, Array< LO > &glCoarseNodesPerDir, Array< GO > &ghostGIDs, Array< GO > &coarseNodesGIDs, Array< GO > &colGIDs, GO &gNumCoarseNodes, LO &lNumCoarseNodes, ArrayRCP< Array< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > > coarseNodes, Array< int > &boundaryFlags, RCP< NodesIDs > ghostedCoarseNodes) const
Class that holds all level-specific information.
Factory that provides an interface for a concrete implementation of a prolongation operator.
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar