MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Zoltan2Interface_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_ZOLTAN2INTERFACE_DEF_HPP
11#define MUELU_ZOLTAN2INTERFACE_DEF_HPP
12
13#include <sstream>
14#include <set>
15
17#if defined(HAVE_MUELU_ZOLTAN2) && defined(HAVE_MPI)
18
19#include <Zoltan2_XpetraMultiVectorAdapter.hpp>
20#include <Zoltan2_XpetraCrsGraphAdapter.hpp>
21#include <Zoltan2_PartitioningProblem.hpp>
22
23#include <Teuchos_Utils.hpp>
24#include <Teuchos_DefaultMpiComm.hpp>
25#include <Teuchos_OpaqueWrapper.hpp>
26
27#include "MueLu_Level.hpp"
28#include "MueLu_Exceptions.hpp"
29#include "MueLu_Monitor.hpp"
30
31namespace MueLu {
32
33template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
35 defaultZoltan2Params = rcp(new ParameterList());
36 defaultZoltan2Params->set("algorithm", "multijagged");
37 defaultZoltan2Params->set("partitioning_approach", "partition");
38
39 // Improve scaling for communication bound algorithms by premigrating
40 // coordinates to a subset of processors.
41 // For more information, see Github issue #1538
42 defaultZoltan2Params->set("mj_premigration_option", 1);
43}
44
45template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
47 RCP<ParameterList> validParamList = rcp(new ParameterList());
48
49 validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Factory of the matrix A");
50 validParamList->set<RCP<const FactoryBase> >("number of partitions", Teuchos::null, "Instance of RepartitionHeuristicFactory.");
51 validParamList->set<RCP<const FactoryBase> >("Coordinates", Teuchos::null, "Factory of the coordinates");
52 validParamList->set<RCP<const ParameterList> >("ParameterList", Teuchos::null, "Zoltan2 parameters");
53 validParamList->set<RCP<const FactoryBase> >("repartition: heuristic target rows per process", Teuchos::null, "Factory for number of rows per process to use with MultiJagged");
54
55 return validParamList;
56}
57
58template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
60 Input(currentLevel, "A");
61 Input(currentLevel, "number of partitions");
62 const ParameterList& pL = GetParameterList();
63 // We do this dance, because we don't want "ParameterList" to be marked as used.
64 // Is there a better way?
65 Teuchos::ParameterEntry entry = pL.getEntry("ParameterList");
66 RCP<const Teuchos::ParameterList> providedList = Teuchos::any_cast<RCP<const Teuchos::ParameterList> >(entry.getAny(false));
67 if (providedList != Teuchos::null && providedList->isType<std::string>("algorithm")) {
68 const std::string algo = providedList->get<std::string>("algorithm");
69 if (algo == "multijagged") {
70 Input(currentLevel, "Coordinates");
71 Input(currentLevel, "repartition: heuristic target rows per process");
72 } else if (algo == "rcb") {
73 Input(currentLevel, "Coordinates");
74 }
75 } else {
76 Input(currentLevel, "repartition: heuristic target rows per process");
77 Input(currentLevel, "Coordinates");
78 }
79}
80
81template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
83 FactoryMonitor m(*this, "Build", level);
84
85 typedef typename Teuchos::ScalarTraits<SC>::coordinateType real_type;
86 typedef typename Xpetra::MultiVector<real_type, LO, GO, NO> RealValuedMultiVector;
88
89 RCP<Matrix> A = Get<RCP<Matrix> >(level, "A");
90
91 int numParts = Get<int>(level, "number of partitions");
92 if (numParts == 1 || numParts == -1) {
93 // Single processor, decomposition is trivial: all zeros
94 RCP<const Map> rowMap = A->getRowMap();
95 RCP<Xpetra::Vector<GO, LO, GO, NO> > decomposition = Xpetra::VectorFactory<GO, LO, GO, NO>::Build(rowMap, true);
96 Set(level, "Partition", decomposition);
97 return;
98 } /* else if (numParts == -1) {
99 // No repartitioning
100 RCP<Xpetra::Vector<GO,LO,GO,NO> > decomposition = Teuchos::null; //Xpetra::VectorFactory<GO, LO, GO, NO>::Build(rowMap, true);
101 //decomposition->putScalar(Teuchos::as<Scalar>(rowMap->getComm()->getRank()));
102 Set(level, "Partition", decomposition);
103 return;
104 }*/
105
106 const ParameterList& pL = GetParameterList();
107
108 RCP<const ParameterList> providedList = pL.get<RCP<const ParameterList> >("ParameterList");
109 ParameterList Zoltan2Params;
110 if (providedList != Teuchos::null)
111 Zoltan2Params = *providedList;
112
113 // Merge defalt Zoltan2 parameters with user provided
114 // If default and user parameters contain the same parameter name, user one is always preferred
115 for (ParameterList::ConstIterator param = defaultZoltan2Params->begin(); param != defaultZoltan2Params->end(); param++) {
116 const std::string& pName = defaultZoltan2Params->name(param);
117 if (!Zoltan2Params.isParameter(pName))
118 Zoltan2Params.setEntry(pName, defaultZoltan2Params->getEntry(pName));
119 }
120 Zoltan2Params.set("num_global_parts", Teuchos::as<int>(numParts));
121
122 const std::string& algo = Zoltan2Params.get<std::string>("algorithm");
123
124 if (algo == "multijagged" && !Zoltan2Params.isParameter("mj_premigration_coordinate_count")) {
125 LO heuristicTargetRowsPerProcess = Get<LO>(level, "repartition: heuristic target rows per process");
126 Zoltan2Params.set("mj_premigration_coordinate_count", heuristicTargetRowsPerProcess);
127 }
128
129 GetOStream(Runtime0) << "Zoltan2 parameters:\n----------\n"
130 << Zoltan2Params << "----------" << std::endl;
131
132 RCP<RealValuedMultiVector> coords;
133
134 if (algo == "multijagged" || algo == "rcb")
135 coords = Get<RCP<RealValuedMultiVector> >(level, "Coordinates");
136
137 const std::string debuggingFile = "mj_debug.lvl_" + std::to_string(level.GetLevelID());
138 RCP<typename Util::GOVector> decomposition;
139 {
140 SubFactoryMonitor m1(*this, "Zoltan2 " + toString(algo), level);
141 decomposition = Util::ComputeDecomposition(numParts, A, coords, Zoltan2Params, debuggingFile);
142 }
143
144 Set(level, "Partition", decomposition);
145}
146
147template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
148RCP<typename Zoltan2Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::GOVector>
149Zoltan2Utilities<Scalar, LocalOrdinal, GlobalOrdinal, Node>::ComputeDecomposition(int numPartitions, RCP<Matrix>& A, RCP<RealValuedMultiVector> coords, Teuchos::ParameterList& Zoltan2Params, std::string debuggingFile) {
150 Zoltan2Params.set("num_global_parts", numPartitions);
151
152 RCP<const Map> rowMap = A->getRowMap();
153 LO blkSize = A->GetFixedBlockSize();
154
155 const std::string& algo = Zoltan2Params.get<std::string>("algorithm");
156
157 if (algo == "multijagged" || algo == "rcb") {
158 TEUCHOS_ASSERT(coords);
159 RCP<const Map> map = coords->getMap();
160 GO numElements = map->getLocalNumElements();
161
162 // Check that the number of local coordinates is consistent with the #rows in A
163 TEUCHOS_TEST_FOR_EXCEPTION(rowMap->getLocalNumElements() / blkSize != coords->getLocalLength(), Exceptions::Incompatible,
164 "Coordinate vector length (" + toString(coords->getLocalLength()) << " is incompatible with number of block rows in A (" + toString(rowMap->getLocalNumElements() / blkSize) + "The vector length should be the same as the number of mesh points.");
165#ifdef HAVE_MUELU_DEBUG
166 GO indexBase = rowMap->getIndexBase();
167 // Make sure that logical blocks in row map coincide with logical nodes in coordinates map
168 ArrayView<const GO> rowElements = rowMap->getLocalElementList();
169 ArrayView<const GO> coordsElements = map->getLocalElementList();
170 for (LO i = 0; i < Teuchos::as<LO>(numElements); i++)
171 TEUCHOS_TEST_FOR_EXCEPTION((coordsElements[i] - indexBase) * blkSize + indexBase != rowElements[i * blkSize],
172 Exceptions::RuntimeError, "i = " << i << ", coords GID = " << coordsElements[i] << ", row GID = " << rowElements[i * blkSize] << ", blkSize = " << blkSize << std::endl);
173#endif
174
175 typedef Zoltan2::XpetraMultiVectorAdapter<RealValuedMultiVector> InputAdapterType;
176 typedef Zoltan2::PartitioningProblem<InputAdapterType> ProblemType;
177
178 Array<real_type> weightsPerRow(numElements);
179 for (LO i = 0; i < numElements; i++) {
180 weightsPerRow[i] = 0.0;
181
182 for (LO j = 0; j < blkSize; j++) {
183 weightsPerRow[i] += A->getNumEntriesInLocalRow(i * blkSize + j);
184 }
185 }
186
187 const bool writeZoltan2DebuggingFiles = Zoltan2Params.get("mj_debug", false);
188 Zoltan2Params.remove("mj_debug");
189
190 std::vector<int> strides;
191 std::vector<const real_type*> weights(1, weightsPerRow.getRawPtr());
192
193 RCP<const Teuchos::MpiComm<int> > dupMpiComm = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(rowMap->getComm()->duplicate());
194 RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > zoltanComm = dupMpiComm->getRawMpiComm();
195
196 InputAdapterType adapter(coords, weights, strides);
197 RCP<ProblemType> problem(new ProblemType(&adapter, &Zoltan2Params, (*zoltanComm)()));
198
199 {
200 if (writeZoltan2DebuggingFiles)
201 adapter.generateFiles(debuggingFile.c_str(), *(rowMap->getComm()));
202 problem->solve();
203 }
204
205 RCP<Xpetra::Vector<GO, LO, GO, NO> > decomposition = Xpetra::VectorFactory<GO, LO, GO, NO>::Build(rowMap, false);
206 ArrayRCP<GO> decompEntries = decomposition->getDataNonConst(0);
207
208 const typename InputAdapterType::part_t* parts = problem->getSolution().getPartListView();
209
210 for (GO i = 0; i < numElements; i++) {
211 int partNum = parts[i];
212
213 for (LO j = 0; j < blkSize; j++)
214 decompEntries[i * blkSize + j] = partNum;
215 }
216
217 return decomposition;
218 }
219
220 GO numElements = rowMap->getLocalNumElements();
221
222 typedef Zoltan2::XpetraCrsGraphAdapter<CrsGraph> InputAdapterType;
223 typedef Zoltan2::PartitioningProblem<InputAdapterType> ProblemType;
224
225 RCP<const Teuchos::MpiComm<int> > dupMpiComm = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(rowMap->getComm()->duplicate());
226 RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > zoltanComm = dupMpiComm->getRawMpiComm();
227
228 InputAdapterType adapter(A->getCrsGraph());
229 RCP<ProblemType> problem(new ProblemType(&adapter, &Zoltan2Params, (*zoltanComm)()));
230
231 problem->solve();
232
233 RCP<Xpetra::Vector<GO, LO, GO, NO> > decomposition = Xpetra::VectorFactory<GO, LO, GO, NO>::Build(rowMap, false);
234 ArrayRCP<GO> decompEntries = decomposition->getDataNonConst(0);
235
236 const typename InputAdapterType::part_t* parts = problem->getSolution().getPartListView();
237
238 // For blkSize > 1, ignore solution for every row but the first ones in a block.
239 for (GO i = 0; i < numElements / blkSize; i++) {
240 int partNum = parts[i * blkSize];
241
242 for (LO j = 0; j < blkSize; j++)
243 decompEntries[i * blkSize + j] = partNum;
244 }
245
246 return decomposition;
247}
248
249} // namespace MueLu
250
251#endif // if defined(HAVE_MUELU_ZOLTAN2) && defined(HAVE_MPI)
252
253#endif // MUELU_ZOLTAN2INTERFACE_DEF_HPP
Exception throws to report incompatible objects (like maps).
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
int GetLevelID() const
Return level number.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
void DeclareInput(Level &currentLevel) 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.
void Build(Level &currentLevel) const
Build an object with this factory.
Namespace for MueLu classes and methods.
@ Runtime0
One-liner description of what is happening.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
static RCP< GOVector > ComputeDecomposition(int numPartitions, RCP< Matrix > &A, RCP< RealValuedMultiVector > coords, Teuchos::ParameterList &Zoltan2Params, std::string debuggingFile="")