46#ifndef MUELU_ZOLTAN2INTERFACE_DECL_HPP
47#define MUELU_ZOLTAN2INTERFACE_DECL_HPP
51#if defined(HAVE_MUELU_ZOLTAN2) && defined(HAVE_MPI)
53#include <Xpetra_Matrix.hpp>
54#include <Xpetra_MultiVectorFactory.hpp>
55#include <Xpetra_VectorFactory.hpp>
64#if defined(HAVE_MUELU_ZOLTAN)
65#include "MueLu_ZoltanInterface.hpp"
119#undef MUELU_ZOLTAN2INTERFACE_SHORT
152#ifdef HAVE_MUELU_EPETRA
154#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
155 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
157#if defined(HAVE_MUELU_ZOLTAN)
165 typedef Xpetra::EpetraNode
Node;
166#undef MUELU_ZOLTAN2INTERFACE_SHORT
169 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType
real_type;
173 level_ = rcp(
new Level());
176 level_->SetLevelID(1);
179 zoltanInterface_ = Teuchos::null;
180 level_ = Teuchos::null;
183 RCP<ParameterList> validParamList = rcp(
new ParameterList());
185 validParamList->set< RCP<const FactoryBase> > (
"A", Teuchos::null,
"Factory of the matrix A");
186 validParamList->set< RCP<const FactoryBase> > (
"Coordinates", Teuchos::null,
"Factory of the coordinates");
187 validParamList->set< RCP<const FactoryBase> > (
"number of partitions", Teuchos::null,
"Instance of RepartitionHeuristicFactory.");
188 validParamList->set< RCP<const ParameterList> > (
"ParameterList", Teuchos::null,
"Zoltan2 parameters");
190 return validParamList;
193 Input(currentLevel,
"A");
194 Input(currentLevel,
"number of partitions");
198 Teuchos::ParameterEntry entry = pL.getEntry(
"ParameterList");
199 RCP<const Teuchos::ParameterList> providedList = Teuchos::any_cast<RCP<const Teuchos::ParameterList> >(entry.getAny(
false));
200 if (providedList != Teuchos::null && providedList->isType<std::string>(
"algorithm")) {
201 const std::string algo = providedList->get<std::string>(
"algorithm");
202 if (algo ==
"multijagged" || algo ==
"rcb")
203 Input(currentLevel,
"Coordinates");
205 Input(currentLevel,
"Coordinates");
208 this->
GetOStream(
Warnings0) <<
"Tpetra does not support <double,int,int,EpetraNode> instantiation, "
209 "switching Zoltan2Interface to ZoltanInterface" << std::endl;
212 level_->Set(
"A",
Get<RCP<Matrix> > (currentLevel,
"A"));
213 level_->Set(
"Coordinates",
Get<RCP<RealValuedMultiVector> >(currentLevel,
"Coordinates"));
214 level_->Set(
"number of partitions", currentLevel.
Get<GO>(
"number of partitions"));
216 level_->Request(
"Partition", zoltanInterface_.get());
217 zoltanInterface_->Build(*level_);
219 RCP<Xpetra::Vector<GO,LO,GO,NO> > decomposition;
220 level_->Get(
"Partition", decomposition, zoltanInterface_.get());
221 Set(currentLevel,
"Partition", decomposition);
237 void Build(Level &level)
const {};
247#define MUELU_ZOLTAN2INTERFACE_SHORT
Kokkos::Compat::KokkosSerialWrapperNode EpetraNode
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Exception throws to report errors in the internal logical of the program.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
void Set(Level &level, const std::string &varName, const T &data) const
Class that holds all level-specific information.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
virtual const Teuchos::ParameterList & GetParameterList() const
Base class for factories that use one level (currentLevel).
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
virtual ~Zoltan2Interface()
RCP< ZoltanInterface > zoltanInterface_
void DeclareInput(Level ¤tLevel) const
Specifies the data that this class needs, and the factories that generate that data.
Teuchos::ScalarTraits< Scalar >::magnitudeType real_type
Xpetra::MultiVector< real_type, LO, GO, NO > RealValuedMultiVector
void Build(Level ¤tLevel) const
Build an object with this factory.
Interface to Zoltan2 library.
void DeclareInput(Level ¤tLevel) const
Specifies the data that this class needs, and the factories that generate that data.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Zoltan2Interface()
Constructor.
void Build(Level ¤tLevel) const
Build an object with this factory.
virtual ~Zoltan2Interface()
Destructor.
RCP< ParameterList > defaultZoltan2Params
Interface to Zoltan library.
Namespace for MueLu classes and methods.
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
@ Warnings0
Important warning messages (one line)