46 #ifndef MUELU_ZOLTAN2INTERFACE_DEF_HPP
47 #define MUELU_ZOLTAN2INTERFACE_DEF_HPP
53 #if defined(HAVE_MUELU_ZOLTAN2) && defined(HAVE_MPI)
55 #include <Zoltan2_XpetraMultiVectorAdapter.hpp>
56 #include <Zoltan2_XpetraCrsGraphAdapter.hpp>
57 #include <Zoltan2_PartitioningProblem.hpp>
61 #include <Teuchos_OpaqueWrapper.hpp>
66 #include "MueLu_Utilities.hpp"
70 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
73 defaultZoltan2Params->set(
"algorithm",
"multijagged");
74 defaultZoltan2Params->set(
"partitioning_approach",
"partition");
79 defaultZoltan2Params->set(
"mj_premigration_option", 1);
82 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
87 validParamList->
set<
RCP<const FactoryBase> > (
"number of partitions", Teuchos::null,
"Instance of RepartitionHeuristicFactory.");
90 validParamList->
set<
RCP<const FactoryBase> > (
"repartition: heuristic target rows per process", Teuchos::null,
"Factory for number of rows per process to use with MultiJagged");
92 return validParamList;
96 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
98 Input(currentLevel,
"A");
99 Input(currentLevel,
"number of partitions");
105 if (providedList != Teuchos::null && providedList->
isType<std::string>(
"algorithm")) {
106 const std::string algo = providedList->
get<std::string>(
"algorithm");
107 if (algo ==
"multijagged") {
108 Input(currentLevel,
"Coordinates");
109 Input(currentLevel,
"repartition: heuristic target rows per process");
110 }
else if (algo ==
"rcb") {
111 Input(currentLevel,
"Coordinates");
114 Input(currentLevel,
"repartition: heuristic target rows per process");
115 Input(currentLevel,
"Coordinates");
119 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
124 typedef typename Xpetra::MultiVector<real_type,LO,GO,NO> RealValuedMultiVector;
128 LO blkSize = A->GetFixedBlockSize();
130 int numParts = Get<int>(level,
"number of partitions");
131 if (numParts == 1 || numParts == -1) {
134 Set(level,
"Partition", decomposition);
148 if (providedList != Teuchos::null)
149 Zoltan2Params = *providedList;
154 const std::string& pName = defaultZoltan2Params->
name(param);
155 if (!Zoltan2Params.isParameter(pName))
156 Zoltan2Params.setEntry(pName, defaultZoltan2Params->getEntry(pName));
158 Zoltan2Params.set(
"num_global_parts", Teuchos::as<int>(numParts));
160 GetOStream(
Runtime0) <<
"Zoltan2 parameters:\n----------\n" << Zoltan2Params <<
"----------" << std::endl;
162 const std::string& algo = Zoltan2Params.
get<std::string>(
"algorithm");
164 if (algo ==
"multijagged" || algo ==
"rcb") {
168 GO numElements = map->getNodeNumElements();
172 "Coordinate vector length (" +
toString(coords->getLocalLength()) <<
" is incompatible with number of block rows in A ("
173 +
toString(rowMap->getNodeNumElements()/blkSize) +
"The vector length should be the same as the number of mesh points.");
174 #ifdef HAVE_MUELU_DEBUG
175 GO indexBase = rowMap->getIndexBase();
176 GetOStream(
Runtime0) <<
"Checking consistence of row and coordinates maps" << std::endl;
180 for (LO i = 0; i < Teuchos::as<LO>(numElements); i++)
183 <<
", row GID = " << rowElements[i*blkSize] <<
", blkSize = " << blkSize << std::endl);
186 typedef Zoltan2::XpetraMultiVectorAdapter<RealValuedMultiVector> InputAdapterType;
187 typedef Zoltan2::PartitioningProblem<InputAdapterType> ProblemType;
190 for (LO i = 0; i < numElements; i++) {
191 weightsPerRow[i] = 0.0;
193 for (LO j = 0; j < blkSize; j++) {
194 weightsPerRow[i] += A->getNumEntriesInLocalRow(i*blkSize+j);
199 if(algo ==
"multijagged" && !Zoltan2Params.isParameter(
"mj_premigration_coordinate_count")) {
200 LO heuristicTargetRowsPerProcess = Get<LO>(level,
"repartition: heuristic target rows per process");
201 Zoltan2Params.set(
"mj_premigration_coordinate_count", heuristicTargetRowsPerProcess);
203 const bool writeZoltan2DebuggingFiles = Zoltan2Params.
get(
"mj_debug",
false);
204 Zoltan2Params.remove(
"mj_debug");
206 std::vector<int> strides;
207 std::vector<const real_type*> weights(1, weightsPerRow.
getRawPtr());
212 InputAdapterType adapter(coords, weights, strides);
213 RCP<ProblemType> problem(
new ProblemType(&adapter, &Zoltan2Params, (*zoltanComm)()));
217 if (writeZoltan2DebuggingFiles)
218 adapter.generateFiles((
"mj_debug.lvl_"+std::to_string(level.GetLevelID())).c_str(), *(rowMap->getComm()));
223 ArrayRCP<GO> decompEntries = decomposition->getDataNonConst(0);
225 const typename InputAdapterType::part_t * parts = problem->getSolution().getPartListView();
227 for (GO i = 0; i < numElements; i++) {
228 int partNum = parts[i];
230 for (LO j = 0; j < blkSize; j++)
231 decompEntries[i*blkSize + j] = partNum;
234 Set(level,
"Partition", decomposition);
238 GO numElements = rowMap->getNodeNumElements();
240 typedef Zoltan2::XpetraCrsGraphAdapter<CrsGraph> InputAdapterType;
241 typedef Zoltan2::PartitioningProblem<InputAdapterType> ProblemType;
246 InputAdapterType adapter(A->getCrsGraph());
247 RCP<ProblemType> problem(
new ProblemType(&adapter, &Zoltan2Params, (*zoltanComm)()));
255 ArrayRCP<GO> decompEntries = decomposition->getDataNonConst(0);
257 const typename InputAdapterType::part_t * parts = problem->getSolution().getPartListView();
260 for (GO i = 0; i < numElements/blkSize; i++) {
261 int partNum = parts[i*blkSize];
263 for (LO j = 0; j < blkSize; j++)
264 decompEntries[i*blkSize + j] = partNum;
267 Set(level,
"Partition", decomposition);
273 #endif //if defined(HAVE_MUELU_ZOLTAN2) && defined(HAVE_MPI)
275 #endif // MUELU_ZOLTAN2INTERFACE_DEF_HPP
const std::string & name() const
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Timer to be used in factories. Similar to Monitor but with additional timers.
Zoltan2Interface()
Constructor.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
One-liner description of what is happening.
void Build(Level ¤tLevel) const
Build an object with this factory.
Exception throws to report incompatible objects (like maps).
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
params_t::ConstIterator ConstIterator
any & getAny(bool activeQry=true)
bool isType(const std::string &name) const
Exception throws to report errors in the internal logical of the program.
ParameterEntry & getEntry(const std::string &name)
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.