MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_IntrepidPCoarsenFactory_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_IPCFACTORY_DECL_HPP
11 #define MUELU_IPCFACTORY_DECL_HPP
12 
13 #include <string>
14 #include <vector>
15 
16 #include "MueLu_ConfigDefs.hpp"
18 
19 #include "MueLu_Level_fwd.hpp"
21 #include "MueLu_PerfUtils_fwd.hpp"
22 #include "MueLu_PFactory.hpp"
23 #include "MueLu_Utilities_fwd.hpp"
24 
25 #include "Intrepid2_Basis.hpp"
26 
27 #include "Kokkos_DynRankView.hpp"
28 
29 #include <Xpetra_Import_fwd.hpp>
30 
31 namespace MueLu {
32 
73 template <class Scalar = DefaultScalar,
76  class Node = DefaultNode>
78 #undef MUELU_INTREPIDPCOARSENFACTORY_SHORT
79 #include "MueLu_UseShortNames.hpp"
80 
81  public:
82  typedef Kokkos::DynRankView<LocalOrdinal, typename Node::device_type> LOFieldContainer;
83  typedef Kokkos::DynRankView<double, typename Node::device_type> SCFieldContainer;
84  typedef Intrepid2::Basis<typename Node::device_type::execution_space, double, double> Basis; // Hardwired on purpose
85 
87 
88 
93 
96 
98 
100 
102 
103 
104  void DeclareInput(Level &fineLevel, Level &coarseLevel) const;
105 
107 
109 
110 
116  void Build(Level &fineLevel, Level &coarseLevel) const;
117 
118  void BuildP(Level &fineLevel, Level &coarseLevel) const; // Build()
119 
121  private:
123 
124  // NOTE: This is hardwired to double on purpose.
126  const std::vector<bool> &hi_nodeIsOwned,
127  const SCFieldContainer &hi_DofCoords,
128  const std::vector<size_t> &lo_node_in_hi,
129  const Basis &lo_Basis,
130  const std::vector<LocalOrdinal> &hi_to_lo_map,
131  const Teuchos::RCP<const Map> &lo_colMap,
132  const Teuchos::RCP<const Map> &lo_domainMap,
133  const Teuchos::RCP<const Map> &hi_map,
134  Teuchos::RCP<Matrix> &P) const;
136 
137  // NOTE: This is hardwired to double on purpose.
139  const std::vector<bool> &hi_nodeIsOwned,
140  const SCFieldContainer &hi_DofCoords,
141  const LOFieldContainer &lo_elemToHiRepresentativeNode,
142  const Basis &lo_basis,
143  const std::vector<LocalOrdinal> &hi_to_lo_map,
144  const Teuchos::RCP<const Map> &lo_colMap,
145  const Teuchos::RCP<const Map> &lo_domainMap,
146  const Teuchos::RCP<const Map> &hi_map,
147  Teuchos::RCP<Matrix> &P) const;
148 
149 }; // class IntrepidPCoarsenFactory
150 
151 /* Utility functions for use with Intrepid */
152 namespace MueLuIntrepid {
153 
154 template <class Scalar, class KokkosExecutionSpace>
156 
157 template <class Scalar, class KokkosDeviceType>
158 void IntrepidGetLoNodeInHi(const Teuchos::RCP<Intrepid2::Basis<typename KokkosDeviceType::execution_space, Scalar, Scalar> > &hi_basis,
159  const Teuchos::RCP<Intrepid2::Basis<typename KokkosDeviceType::execution_space, Scalar, Scalar> > &lo_basis,
160  std::vector<size_t> &lo_node_in_hi,
161  Kokkos::DynRankView<Scalar, KokkosDeviceType> &hi_DofCoords);
162 
163 template <class LocalOrdinal, class GlobalOrdinal, class Node, class LOFieldContainer>
164 void GenerateLoNodeInHiViaGIDs(const std::vector<std::vector<size_t> > &candidates, const LOFieldContainer &hi_elemToNode,
166  LOFieldContainer &lo_elemToHiRepresentativeNode);
167 
168 template <class LocalOrdinal, class LOFieldContainer>
169 void BuildLoElemToNode(const LOFieldContainer &hi_elemToNode,
170  const std::vector<bool> &hi_nodeIsOwned,
171  const std::vector<size_t> &lo_node_in_hi,
172  const Teuchos::ArrayRCP<const int> &hi_isDirichlet,
173  LOFieldContainer &lo_elemToNode,
174  std::vector<bool> &lo_nodeIsOwned,
175  std::vector<LocalOrdinal> &hi_to_lo_map,
176  int &lo_numOwnedNodes);
177 
178 template <class LocalOrdinal, class LOFieldContainer>
179 void BuildLoElemToNodeViaRepresentatives(const LOFieldContainer &hi_elemToNode,
180  const std::vector<bool> &hi_nodeIsOwned,
181  const LOFieldContainer &lo_elemToHiRepresentativeNode,
182  LOFieldContainer &lo_elemToNode,
183  std::vector<bool> &lo_nodeIsOwned,
184  std::vector<LocalOrdinal> &hi_to_lo_map,
185  int &lo_numOwnedNodes);
186 
187 template <class LocalOrdinal, class GlobalOrdinal, class Node>
188 void GenerateColMapFromImport(const Xpetra::Import<LocalOrdinal, GlobalOrdinal, Node> &hi_importer, const std::vector<LocalOrdinal> &hi_to_lo_map, const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> &lo_domainMap, const size_t &lo_columnMapLength, RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > &lo_columnMap);
189 
190 template <class Basis, class SCFieldContainer>
191 void GenerateRepresentativeBasisNodes(const Basis &basis, const SCFieldContainer &ReferenceNodeLocations, const double threshold, std::vector<std::vector<size_t> > &representative_node_candidates);
192 
193 // ! Given an element to (global) node map and a basis, determine one global ordinal per geometric entity (vertex, edge, face,
194 // ! interior). On exit, seeds container is of dimension (spaceDim+1), and contains a sorted vector of local ordinals
195 // ! belonging to entities of that dimension. Only locally-owned degrees of freedom (as determined by rowMap and columnMap)
196 // ! will be stored in seeds.
197 template <class Basis, class LOFieldContainer, class LocalOrdinal, class GlobalOrdinal, class Node>
198 void FindGeometricSeedOrdinals(Teuchos::RCP<Basis> basis, const LOFieldContainer &elementToNodeMap,
199  std::vector<std::vector<LocalOrdinal> > &seeds,
202 
203 } // namespace MueLuIntrepid
204 } // namespace MueLu
205 
206 #define MUELU_INTREPIDPCOARSENFACTORY_SHORT
207 #endif // MUELU_INTREPIDPCOARSENFACTORY_DECL_HPP
void GenerateLoNodeInHiViaGIDs(const std::vector< std::vector< size_t > > &candidates, const LOFieldContainer &hi_elemToNode, RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &hi_columnMap, LOFieldContainer &lo_elemToHiRepresentativeNode)
MueLu::DefaultLocalOrdinal LocalOrdinal
Teuchos::RCP< Intrepid2::Basis< KokkosExecutionSpace, Scalar, Scalar > > BasisFactory(const std::string &name, int &degree)
void FindGeometricSeedOrdinals(Teuchos::RCP< Basis > basis, const LOFieldContainer &elementToNodeMap, std::vector< std::vector< LocalOrdinal > > &seeds, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &rowMap, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &columnMap)
void IntrepidGetLoNodeInHi(const Teuchos::RCP< Intrepid2::Basis< typename KokkosDeviceType::execution_space, Scalar, Scalar > > &hi_basis, const Teuchos::RCP< Intrepid2::Basis< typename KokkosDeviceType::execution_space, Scalar, Scalar > > &lo_basis, std::vector< size_t > &lo_node_in_hi, Kokkos::DynRankView< Scalar, KokkosDeviceType > &hi_DofCoords)
Factory for building transfer operators based on coarsening in polynomial degree, following the Intre...
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
MueLu::DefaultNode Node
Kokkos::DynRankView< LocalOrdinal, typename Node::device_type > LOFieldContainer
IntrepidPCoarsenFactory()
Constructor. User can supply a factory for generating the tentative prolongator.
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 BuildLoElemToNode(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const std::vector< size_t > &lo_node_in_hi, const Teuchos::ArrayRCP< const int > &hi_isDirichlet, LOFieldContainer &lo_elemToNode, std::vector< bool > &lo_nodeIsOwned, std::vector< LocalOrdinal > &hi_to_lo_map, int &lo_numOwnedNodes)
Kokkos::DynRankView< double, typename Node::device_type > SCFieldContainer
void BuildLoElemToNodeViaRepresentatives(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const LOFieldContainer &lo_elemToHiRepresentativeNode, LOFieldContainer &lo_elemToNode, std::vector< bool > &lo_nodeIsOwned, std::vector< LocalOrdinal > &hi_to_lo_map, int &lo_numOwnedNodes)
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void GenerateLinearCoarsening_pn_kirby_to_p1(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const SCFieldContainer &hi_DofCoords, const std::vector< size_t > &lo_node_in_hi, const Basis &lo_Basis, const std::vector< LocalOrdinal > &hi_to_lo_map, const Teuchos::RCP< const Map > &lo_colMap, const Teuchos::RCP< const Map > &lo_domainMap, const Teuchos::RCP< const Map > &hi_map, Teuchos::RCP< Matrix > &P) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void GenerateColMapFromImport(const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &hi_importer, const std::vector< LocalOrdinal > &hi_to_lo_map, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &lo_domainMap, const size_t &lo_columnMapLength, RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &lo_columnMap)
Intrepid2::Basis< typename Node::device_type::execution_space, double, double > Basis
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
void GenerateRepresentativeBasisNodes(const Basis &basis, const SCFieldContainer &ReferenceNodeLocations, const double threshold, std::vector< std::vector< size_t > > &representative_node_candidates)
void GenerateLinearCoarsening_pn_kirby_to_pm(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const SCFieldContainer &hi_DofCoords, const LOFieldContainer &lo_elemToHiRepresentativeNode, const Basis &lo_basis, const std::vector< LocalOrdinal > &hi_to_lo_map, const Teuchos::RCP< const Map > &lo_colMap, const Teuchos::RCP< const Map > &lo_domainMap, const Teuchos::RCP< const Map > &hi_map, Teuchos::RCP< Matrix > &P) const
Factory that provides an interface for a concrete implementation of a prolongation operator...