46 #ifndef MUELU_BLACKBOXPFACTORY_DECL_HPP 47 #define MUELU_BLACKBOXPFACTORY_DECL_HPP 49 #include <Teuchos_SerialDenseVector.hpp> 51 #include <Xpetra_MultiVector.hpp> 52 #include <Xpetra_Matrix_fwd.hpp> 55 #include "MueLu_PFactory.hpp" 62 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
116 template <
class Scalar =
double,
class LocalOrdinal =
int,
class GlobalOrdinal = LocalOrdinal,
class Node = KokkosClassic::DefaultNode::DefaultNodeType>
118 #undef MUELU_BLACKBOXPFACTORY_SHORT 173 void GetGeometricData(RCP<Xpetra::MultiVector<double,LO,GO,NO> >& coordinates,
174 const Array<LO> coarseRate,
const Array<GO> gFineNodesPerDir,
175 const Array<LO> lFineNodesPerDir,
const LO BlkSize, Array<GO>& gIndices,
176 Array<LO>& myOffset, Array<bool>& ghostInterface, Array<LO>& endRate,
177 Array<GO>& gCoarseNodesPerDir, Array<LO>& lCoarseNodesPerDir,
178 Array<LO>& glCoarseNodesPerDir, Array<GO>& ghostGIDs,
179 Array<GO>& coarseNodesGIDs, Array<GO>& colGIDs, GO& gNumCoarseNodes,
180 LO& lNumCoarseNodes, ArrayRCP<Array<double> > coarseNodes,
181 Array<int>& boundaryFlags, RCP<NodesIDs> ghostedCoarseNodes)
const;
184 const Array<LO> endRate,
const LO BlkSize,
const Array<LO> elemInds,
185 const Array<LO> lCoarseElementsPerDir,
186 const LO numDimensions,
const Array<LO> lFineNodesPerDir,
187 const Array<GO> gFineNodesPerDir,
const Array<GO> gIndices,
188 const Array<LO> lCoarseNodesPerDir,
const Array<bool> ghostInterface,
189 const Array<int> elementFlags,
const std::string stencilType,
190 const std::string blockStrategy,
const Array<LO> elementNodesPerDir,
191 const LO numNodesInElement,
const Array<GO> colGIDs,
192 Teuchos::SerialDenseMatrix<LO,SC>& Pi,
193 Teuchos::SerialDenseMatrix<LO,SC>& Pf,
194 Teuchos::SerialDenseMatrix<LO,SC>& Pe,
195 Array<LO>& dofType, Array<LO>& lDofInd)
const;
197 void CollapseStencil(
const int type,
const int orientation,
const int collapseFlags[3],
198 Array<SC>& stencil)
const ;
200 void FormatStencil(
const LO BlkSize,
const Array<bool> ghostInterface,
const LO ie,
201 const LO je,
const LO ke,
const ArrayView<const SC> rowValues,
202 const Array<LO> elementNodesPerDir,
const int collapseFlags[3],
203 const std::string stencilType, Array<SC>& stencil)
const;
205 void GetNodeInfo(
const LO ie,
const LO je,
const LO ke,
const Array<LO> elementNodesPerDir,
206 int* type, LO& ind,
int* orientation)
const;
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;
218 #define MUELU_BLACKBOXPFACTORY_SHORT 219 #endif // MUELU_BLACKBOXPFACTORY_DECL_HPP virtual ~BlackBoxPFactory()
Destructor.
BlackBoxPFactory()
Constructor.
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 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 GetGeometricData(RCP< Xpetra::MultiVector< double, 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< double > > coarseNodes, Array< int > &boundaryFlags, RCP< NodesIDs > ghostedCoarseNodes) const
Namespace for MueLu classes and methods.
Prolongator factory performing geometric coarsening.
Class that holds all level-specific information.
void GetNodeInfo(const LO ie, const LO je, const LO ke, const Array< LO > elementNodesPerDir, int *type, LO &ind, int *orientation) const
void CollapseStencil(const int type, const int orientation, const int collapseFlags[3], Array< SC > &stencil) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
std::vector< GO > colInds
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
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
Factory that provides an interface for a concrete implementation of a prolongation operator...
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.