10#ifndef MUELU_IPCFACTORY_DECL_HPP
11#define MUELU_IPCFACTORY_DECL_HPP
22#include "MueLu_PFactory.hpp"
25#include "Intrepid2_Basis.hpp"
27#include "Kokkos_DynRankView.hpp"
29#include <Xpetra_Import_fwd.hpp>
78#undef MUELU_INTREPIDPCOARSENFACTORY_SHORT
82 typedef Kokkos::DynRankView<LocalOrdinal, typename Node::device_type>
LOFieldContainer;
84 typedef Intrepid2::Basis<typename Node::device_type::execution_space, double, double>
Basis;
126 const std::vector<bool> &hi_nodeIsOwned,
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;
139 const std::vector<bool> &hi_nodeIsOwned,
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;
152namespace MueLuIntrepid {
154template <
class Scalar,
class KokkosExecutionSpace>
155Teuchos::RCP<Intrepid2::Basis<KokkosExecutionSpace, Scalar, Scalar> >
BasisFactory(
const std::string &name,
int °ree);
157template <
class Scalar,
class KokkosDeviceType>
158void 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);
163template <
class LocalOrdinal,
class GlobalOrdinal,
class Node,
class LOFieldContainer>
164void GenerateLoNodeInHiViaGIDs(
const std::vector<std::vector<size_t> > &candidates,
const LOFieldContainer &hi_elemToNode,
165 RCP<
const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > &hi_columnMap,
166 LOFieldContainer &lo_elemToHiRepresentativeNode);
168template <
class LocalOrdinal,
class LOFieldContainer>
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);
178template <
class LocalOrdinal,
class LOFieldContainer>
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);
187template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
188void 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);
190template <
class Basis,
class SCFieldContainer>
191void GenerateRepresentativeBasisNodes(
const Basis &basis,
const SCFieldContainer &ReferenceNodeLocations,
const double threshold, std::vector<std::vector<size_t> > &representative_node_candidates);
197template <
class Basis,
class LOFieldContainer,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
199 std::vector<std::vector<LocalOrdinal> > &seeds,
200 const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> &rowMap,
201 const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> &columnMap);
206#define MUELU_INTREPIDPCOARSENFACTORY_SHORT
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Factory for building transfer operators based on coarsening in polynomial degree, following the Intre...
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Kokkos::DynRankView< LocalOrdinal, typename Node::device_type > LOFieldContainer
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
IntrepidPCoarsenFactory()
Constructor. User can supply a factory for generating the tentative prolongator.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
Kokkos::DynRankView< double, typename Node::device_type > SCFieldContainer
virtual ~IntrepidPCoarsenFactory()
Destructor.
Intrepid2::Basis< typename Node::device_type::execution_space, double, double > Basis
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
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
Class that holds all level-specific information.
Factory that provides an interface for a concrete implementation of a prolongation operator.
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)
void GenerateRepresentativeBasisNodes(const Basis &basis, const SCFieldContainer &ReferenceNodeLocations, const double threshold, std::vector< std::vector< size_t > > &representative_node_candidates)
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 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)
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)
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)
Teuchos::RCP< Intrepid2::Basis< KokkosExecutionSpace, Scalar, Scalar > > BasisFactory(const std::string &name, int °ree)
Namespace for MueLu classes and methods.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar