46 #ifndef MUELU_IPCFACTORY_DECL_HPP
47 #define MUELU_IPCFACTORY_DECL_HPP
58 #include "MueLu_PFactory.hpp"
63 #include "Intrepid2_Basis.hpp"
67 #include <Xpetra_Import.hpp>
116 #undef MUELU_INTREPIDPCOARSENFACTORY_SHORT
122 typedef Intrepid2::Basis<typename Node::device_type::execution_space,double,double>
Basis;
164 const std::vector<bool> & hi_nodeIsOwned,
166 const std::vector<size_t> &lo_node_in_hi,
167 const Basis &lo_Basis,
168 const std::vector<LocalOrdinal> & hi_to_lo_map,
177 const std::vector<bool> & hi_nodeIsOwned,
180 const Basis &lo_basis,
181 const std::vector<LocalOrdinal> & hi_to_lo_map,
192 namespace MueLuIntrepid {
194 template<
class Scalar,
class KokkosExecutionSpace>
197 template<
class Scalar,
class KokkosDeviceType>
199 const Teuchos::RCP<Intrepid2::Basis<typename KokkosDeviceType::execution_space,Scalar,Scalar> > &lo_basis,
200 std::vector<size_t> & lo_node_in_hi,
201 Kokkos::DynRankView<Scalar,KokkosDeviceType> & hi_DofCoords);
203 template<
class LocalOrdinal,
class GlobalOrdinal,
class Node,
class LOFieldContainer>
204 void GenerateLoNodeInHiViaGIDs(
const std::vector<std::vector<size_t> > & candidates,
const LOFieldContainer & hi_elemToNode,
205 RCP<
const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > & hi_columnMap,
206 LOFieldContainer & lo_elemToHiRepresentativeNode);
208 template <
class LocalOrdinal,
class LOFieldContainer>
210 const std::vector<bool> & hi_nodeIsOwned,
211 const std::vector<size_t> & lo_node_in_hi,
213 LOFieldContainer & lo_elemToNode,
214 std::vector<bool> & lo_nodeIsOwned,
215 std::vector<LocalOrdinal> & hi_to_lo_map,
216 int & lo_numOwnedNodes);
218 template <
class LocalOrdinal,
class LOFieldContainer>
220 const std::vector<bool> & hi_nodeIsOwned,
221 const LOFieldContainer & lo_elemToHiRepresentativeNode,
222 LOFieldContainer & lo_elemToNode,
223 std::vector<bool> & lo_nodeIsOwned,
224 std::vector<LocalOrdinal> & hi_to_lo_map,
225 int & lo_numOwnedNodes);
228 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
229 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);
232 template<
class Basis,
class SCFieldContainer>
233 void GenerateRepresentativeBasisNodes(
const Basis & basis,
const SCFieldContainer & ReferenceNodeLocations,
const double threshold, std::vector<std::vector<size_t> > & representative_node_candidates);
239 template<
class Basis,
class LOFieldContainer,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
241 std::vector<std::vector<LocalOrdinal> > &seeds,
242 const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> &rowMap,
243 const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> &columnMap);
248 #define MUELU_INTREPIDPCOARSENFACTORY_SHORT
249 #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 °ree)
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)
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
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...
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
IntrepidPCoarsenFactory()
Constructor. User can supply a factory for generating the tentative prolongator.
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
Kokkos::DynRankView< double, typename Node::device_type > SCFieldContainer
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)
virtual ~IntrepidPCoarsenFactory()
Destructor.
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.
Kokkos::DynRankView< LocalOrdinal, typename Node::device_type > LOFieldContainer
Intrepid2::Basis< typename Node::device_type::execution_space, double, double > Basis
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)
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...