MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_GeneralGeometricPFactory_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_GENERALGEOMETRICPFACTORY_DECL_HPP
11 #define MUELU_GENERALGEOMETRICPFACTORY_DECL_HPP
12 
14 
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 
24 namespace MueLuTests {
25 // Forward declaration of friend tester class used to UnitTest GeneralGeometricPFactory
26 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
28 } // namespace MueLuTests
29 
30 namespace MueLu {
31 
79 template <class Scalar = DefaultScalar,
82  class Node = DefaultNode>
84 #undef MUELU_GENERALGEOMETRICPFACTORY_SHORT
85 #include "MueLu_UseShortNames.hpp"
86 
87  public:
89 
91 
92 
95 
99 
101 
103 
104 
105  void DeclareInput(Level& fineLevel, Level& coarseLevel) const;
106 
108 
110 
111 
112  void Build(Level& fineLevel, Level& coarseLevel) const;
113  void BuildP(Level& fineLevel, Level& coarseLevel) const;
114 
116 
117  private:
118  struct GeometricData {
119  // Geometric algorithm require a copious amount of data to be passed around so this struct
120  // will reduce the amount of input/output parameters of methods in the class. Additionally
121  // the struct can be rewritten to accomodate constraints of Kokkos/CUDA data types
122 
123  std::string meshLayout = "Global Lexicographic";
131  std::vector<std::vector<GO> > meshData; // These are sorted later so they are in std::vector
132  bool ghostInterface[6] = {false};
133 
135  coarseRate.resize(3);
136  endRate.resize(3);
139  offsets.resize(6);
141  startIndices.resize(6);
145  }
146  };
147 
148  struct NodesIDs {
149  // This small struct just carries basic data associated with coarse nodes that is needed
150  // to compute colMapP and to fillComplete P,
151 
155  std::vector<GO> colInds;
156  };
157 
158  struct NodeID {
159  // This small struct is similar to the one above but only for one node.
160  // It is used to create a vector of NodeID that can easily be sorted
161 
163  int PID;
165  };
166 
167  void MeshLayoutInterface(const int interpolationOrder, const LO blkSize,
168  RCP<const Map> fineCoordsMap, RCP<GeometricData> myGeometry,
169  RCP<NodesIDs> ghostedCoarseNodes,
170  Array<Array<GO> >& lCoarseNodesGIDs) const;
171 
172  void GetCoarsePoints(const int interpolationOrder, const LO blkSize,
173  RCP<const Map> fineCoordsMap, RCP<GeometricData> myGeometry,
174  RCP<NodesIDs> ghostedCoarseNodes,
175  Array<Array<GO> >& lCoarseNodesGIDs) const;
176 
179  const LO nnzP, const LO dofsPerNode,
180  RCP<const Map>& stridedDomainMapP,
181  RCP<Matrix>& Amat, RCP<Matrix>& P,
183  RCP<NodesIDs> ghostedCoarseNodes, Array<Array<GO> > coarseNodesGIDs,
184  int interpolationOrder) const;
185 
186  void ComputeStencil(const LO numDimension, const Array<GO> currentNodeIndices,
187  const Array<GO> coarseNodeIndices, const LO rate[3],
188  const Array<Array<typename Teuchos::ScalarTraits<Scalar>::coordinateType> > coord, const int interpolationOrder,
189  std::vector<double>& stencil) const;
190 
191  void ComputeConstantInterpolationStencil(const LO numDimension,
192  const Array<GO> currentNodeIndices,
193  const Array<GO> coarseNodeIndices,
194  const LO rate[3], std::vector<double>& stencil) const;
195 
196  void ComputeLinearInterpolationStencil(const LO numDimension, const Array<Array<typename Teuchos::ScalarTraits<Scalar>::coordinateType> > coord,
197  std::vector<double>& stencil) const;
198  void GetInterpolationFunctions(const LO numDimension,
200  double functions[4][8]) const;
201 
202  void sh_sort_permute(
203  const typename Teuchos::Array<LocalOrdinal>::iterator& first1,
204  const typename Teuchos::Array<LocalOrdinal>::iterator& last1,
205  const typename Teuchos::Array<LocalOrdinal>::iterator& first2,
206  const typename Teuchos::Array<LocalOrdinal>::iterator& last2) const;
207 
208  void sh_sort2(
209  const typename Teuchos::Array<LocalOrdinal>::iterator& first1,
210  const typename Teuchos::Array<LocalOrdinal>::iterator& last1,
211  const typename Teuchos::Array<LocalOrdinal>::iterator& first2,
212  const typename Teuchos::Array<LocalOrdinal>::iterator& last2) const;
213 
214  void GetGIDLocalLexicographic(const GO i, const GO j, const GO k,
215  const Array<LO> coarseNodeFineIndices,
216  const RCP<GeometricData> myGeo, const LO myRankIndex, const LO pi,
217  const LO pj, const LO pk,
218  const typename std::vector<std::vector<GO> >::iterator blockStart,
219  const typename std::vector<std::vector<GO> >::iterator blockEnd,
220  GO& myGID, LO& myPID, LO& myLID) const;
221 
222 }; // class GeneralGeometricPFactory
223 
224 } // namespace MueLu
225 
226 #define MUELU_GENERALGEOMETRICPFACTORY_SHORT
227 #endif // MUELU_GENERALGEOMETRICPFACTORY_DECL_HPP
void ComputeLinearInterpolationStencil(const LO numDimension, const Array< Array< typename Teuchos::ScalarTraits< Scalar >::coordinateType > > coord, std::vector< double > &stencil) const
MueLu::DefaultLocalOrdinal LocalOrdinal
void GetCoarsePoints(const int interpolationOrder, const LO blkSize, RCP< const Map > fineCoordsMap, RCP< GeometricData > myGeometry, RCP< NodesIDs > ghostedCoarseNodes, Array< Array< GO > > &lCoarseNodesGIDs) const
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
void MakeGeneralGeometricP(RCP< GeometricData > myGeo, const RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::coordinateType, LO, GO, NO > > &fCoords, const LO nnzP, const LO dofsPerNode, RCP< const Map > &stridedDomainMapP, RCP< Matrix > &Amat, RCP< Matrix > &P, RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::coordinateType, LO, GO, NO > > &cCoords, RCP< NodesIDs > ghostedCoarseNodes, Array< Array< GO > > coarseNodesGIDs, int interpolationOrder) const
GlobalOrdinal GO
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
LocalOrdinal LO
MueLu::DefaultNode Node
Prolongator factory performing geometric coarsening.
void MeshLayoutInterface(const int interpolationOrder, const LO blkSize, RCP< const Map > fineCoordsMap, RCP< GeometricData > myGeometry, RCP< NodesIDs > ghostedCoarseNodes, Array< Array< GO > > &lCoarseNodesGIDs) const
MueLu::DefaultScalar Scalar
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
void sh_sort2(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 ComputeStencil(const LO numDimension, const Array< GO > currentNodeIndices, const Array< GO > coarseNodeIndices, const LO rate[3], const Array< Array< typename Teuchos::ScalarTraits< Scalar >::coordinateType > > coord, const int interpolationOrder, std::vector< double > &stencil) const
void GetGIDLocalLexicographic(const GO i, const GO j, const GO k, const Array< LO > coarseNodeFineIndices, const RCP< GeometricData > myGeo, const LO myRankIndex, const LO pi, const LO pj, const LO pk, const typename std::vector< std::vector< GO > >::iterator blockStart, const typename std::vector< std::vector< GO > >::iterator blockEnd, GO &myGID, LO &myPID, LO &myLID) const
void resize(size_type new_size, const value_type &x=value_type())
void GetInterpolationFunctions(const LO numDimension, const Teuchos::SerialDenseVector< LO, double > parameters, double functions[4][8]) const
void ComputeConstantInterpolationStencil(const LO numDimension, const Array< GO > currentNodeIndices, const Array< GO > coarseNodeIndices, const LO rate[3], std::vector< double > &stencil) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
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.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Node NO
Factory that provides an interface for a concrete implementation of a prolongation operator...
std::vector< T >::iterator iterator