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)
166 #undef MUELU_ZOLTAN2INTERFACE_SHORT
176 level_->SetLevelID(1);
179 zoltanInterface_ = Teuchos::null;
180 level_ = Teuchos::null;
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.");
190 return validParamList;
193 Input(currentLevel,
"A");
194 Input(currentLevel,
"number of partitions");
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"));
214 level_->Set(
"number of partitions", currentLevel.
Get<GO>(
"number of partitions"));
216 level_->Request(
"Partition", zoltanInterface_.get());
217 zoltanInterface_->Build(*level_);
220 level_->Get(
"Partition", decomposition, zoltanInterface_.
get());
221 Set(currentLevel,
"Partition", decomposition);
237 void Build(Level &level)
const {};
239 #endif // HAVE_MUELU_ZOLTAN
243 #endif // HAVE_MUELU_EPETRA
247 #define MUELU_ZOLTAN2INTERFACE_SHORT
248 #endif //if defined(HAVE_MUELU_ZOLTAN2) && defined(HAVE_MPI)
250 #endif // MUELU_ZOLTAN2INTERFACE_DECL_HPP
Important warning messages (one line)
virtual const Teuchos::ParameterList & GetParameterList() const
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
MueLu::DefaultLocalOrdinal LocalOrdinal
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access). Usage: Level->Get< RCP<Matrix> >("A", factory) if factory == NULL => use default factory.
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
RCP< ParameterList > defaultZoltan2Params
T & get(const std::string &name, T def_value)
Zoltan2Interface()
Constructor.
void Build(Level ¤tLevel) const
Build an object with this factory.
Interface to Zoltan library.This interface provides access to partitioning methods in Zoltan...
virtual ~Zoltan2Interface()
Destructor.
void Build(Level ¤tLevel) const
Build an object with this factory.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
Xpetra::MultiVector< real_type, LO, GO, NO > RealValuedMultiVector
Interface to Zoltan2 library.This interface provides access to partitioning methods in Zoltan2...
RCP< ZoltanInterface > zoltanInterface_
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
any & getAny(bool activeQry=true)
bool isType(const std::string &name) const
virtual ~Zoltan2Interface()
Exception throws to report errors in the internal logical of the program.
void DeclareInput(Level ¤tLevel) const
Specifies the data that this class needs, and the factories that generate that data.
Kokkos::Compat::KokkosSerialWrapperNode EpetraNode
ParameterEntry & getEntry(const std::string &name)
Teuchos::ScalarTraits< Scalar >::magnitudeType real_type
Base class for factories that use one level (currentLevel).
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void DeclareInput(Level ¤tLevel) const
Specifies the data that this class needs, and the factories that generate that data.