MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_RebalanceTransferFactory_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
47 #define MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
48 
49 #include <sstream>
50 #include <Teuchos_Tuple.hpp>
51 
52 #include "Xpetra_MultiVector.hpp"
53 #include "Xpetra_MultiVectorFactory.hpp"
54 #include "Xpetra_Vector.hpp"
55 #include "Xpetra_VectorFactory.hpp"
56 #include <Xpetra_Matrix.hpp>
57 #include <Xpetra_MapFactory.hpp>
58 #include <Xpetra_MatrixFactory.hpp>
59 #include <Xpetra_Import.hpp>
60 #include <Xpetra_ImportFactory.hpp>
61 #include <Xpetra_IO.hpp>
62 
64 
65 #include "MueLu_Level.hpp"
66 #include "MueLu_MasterList.hpp"
67 #include "MueLu_Monitor.hpp"
68 #include "MueLu_PerfUtils.hpp"
69 
70 namespace MueLu {
71 
72  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
74  RCP<ParameterList> validParamList = rcp(new ParameterList());
75 
76 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
77  SET_VALID_ENTRY("repartition: rebalance P and R");
78  SET_VALID_ENTRY("repartition: rebalance Nullspace");
79  SET_VALID_ENTRY("transpose: use implicit");
80  SET_VALID_ENTRY("repartition: use subcommunicators");
81 #undef SET_VALID_ENTRY
82 
83  {
85  RCP<validatorType> typeValidator = rcp (new validatorType(Teuchos::tuple<std::string>("Interpolation", "Restriction"), "type"));
86  validParamList->set("type", "Interpolation", "Type of the transfer operator that need to be rebalanced (Interpolation or Restriction)", typeValidator);
87  }
88 
89  validParamList->set< RCP<const FactoryBase> >("P", null, "Factory of the prolongation operator that need to be rebalanced (only used if type=Interpolation)");
90  validParamList->set< RCP<const FactoryBase> >("R", null, "Factory of the restriction operator that need to be rebalanced (only used if type=Restriction)");
91  validParamList->set< RCP<const FactoryBase> >("Nullspace", null, "Factory of the nullspace that need to be rebalanced (only used if type=Interpolation)");
92  validParamList->set< RCP<const FactoryBase> >("Coordinates", null, "Factory of the coordinates that need to be rebalanced (only used if type=Interpolation)");
93  validParamList->set< RCP<const FactoryBase> >("Importer", null, "Factory of the importer object used for the rebalancing");
94  validParamList->set< int > ("write start", -1, "First level at which coordinates should be written to file");
95  validParamList->set< int > ("write end", -1, "Last level at which coordinates should be written to file");
96 
97  // TODO validation: "P" parameter valid only for type="Interpolation" and "R" valid only for type="Restriction". Like so:
98  // if (paramList.isEntry("type") && paramList.get("type) == "Interpolation) {
99  // validParamList->set< RCP<const FactoryBase> >("P", Teuchos::null, "Factory of the prolongation operator that need to be rebalanced (only used if type=Interpolation)");
100 
101  return validParamList;
102  }
103 
104  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
106  const ParameterList& pL = GetParameterList();
107 
108  if (pL.get<std::string>("type") == "Interpolation") {
109  Input(coarseLevel, "P");
110  if (pL.get<bool>("repartition: rebalance Nullspace"))
111  Input(coarseLevel, "Nullspace");
112  if (pL.get< RCP<const FactoryBase> >("Coordinates") != Teuchos::null)
113  Input(coarseLevel, "Coordinates");
114 
115  } else {
116  if (pL.get<bool>("transpose: use implicit") == false)
117  Input(coarseLevel, "R");
118  }
119 
120  Input(coarseLevel, "Importer");
121  }
122 
123  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
125  FactoryMonitor m(*this, "Build", coarseLevel);
126  typedef Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO> xdMV;
127 
128  const ParameterList& pL = GetParameterList();
129 
130  int implicit = !pL.get<bool>("repartition: rebalance P and R");
131  int writeStart = pL.get<int> ("write start");
132  int writeEnd = pL.get<int> ("write end");
133 
134  if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd && IsAvailable(fineLevel, "Coordinates")) {
135  std::string fileName = "coordinates_level_0.m";
136  RCP<xdMV> fineCoords = fineLevel.Get< RCP<xdMV> >("Coordinates");
137  if (fineCoords != Teuchos::null)
138  Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO>::Write(fileName, *fineCoords);
139  }
140 
141  RCP<const Import> importer = Get<RCP<const Import> >(coarseLevel, "Importer");
142  if (implicit) {
143  // Save the importer, we'll need it for solve
144  coarseLevel.Set("Importer", importer, NoFactory::get());
145  }
146 
147  RCP<ParameterList> params = rcp(new ParameterList());
148  if (IsPrint(Statistics2)) {
149  params->set("printLoadBalancingInfo", true);
150  params->set("printCommInfo", true);
151  }
152 
153  std::string transferType = pL.get<std::string>("type");
154  if (transferType == "Interpolation") {
155  RCP<Matrix> originalP = Get< RCP<Matrix> >(coarseLevel, "P");
156 
157  {
158  // This line must be after the Get call
159  SubFactoryMonitor m1(*this, "Rebalancing prolongator", coarseLevel);
160 
161  if (implicit || importer.is_null()) {
162  GetOStream(Runtime0) << "Using original prolongator" << std::endl;
163  Set(coarseLevel, "P", originalP);
164 
165  } else {
166  // P is the transfer operator from the coarse grid to the fine grid.
167  // P must transfer the data from the newly reordered coarse A to the
168  // (unchanged) fine A. This means that the domain map (coarse) of P
169  // must be changed according to the new partition. The range map
170  // (fine) is kept unchanged.
171  //
172  // The domain map of P must match the range map of R. See also note
173  // below about domain/range map of R and its implications for P.
174  //
175  // To change the domain map of P, P needs to be fillCompleted again
176  // with the new domain map. To achieve this, P is copied into a new
177  // matrix that is not fill-completed. The doImport() operation is
178  // just used here to make a copy of P: the importer is trivial and
179  // there is no data movement involved. The reordering actually
180  // happens during the fillComplete() with domainMap == importer->getTargetMap().
181  RCP<Matrix> rebalancedP = originalP;
182  RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(originalP);
183  TEUCHOS_TEST_FOR_EXCEPTION(crsOp == Teuchos::null, Exceptions::BadCast, "Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
184 
185  RCP<CrsMatrix> rebalancedP2 = crsOp->getCrsMatrix();
186  TEUCHOS_TEST_FOR_EXCEPTION(rebalancedP2 == Teuchos::null, std::runtime_error, "Xpetra::CrsMatrixWrap doesn't have a CrsMatrix");
187 
188  {
189  SubFactoryMonitor subM(*this, "Rebalancing prolongator -- fast map replacement", coarseLevel);
190 
191  RCP<const Import> newImporter;
192  {
193  SubFactoryMonitor(*this, "Import construction", coarseLevel);
194  newImporter = ImportFactory::Build(importer->getTargetMap(), rebalancedP->getColMap());
195  }
196  rebalancedP2->replaceDomainMapAndImporter(importer->getTargetMap(), newImporter);
197  }
198 
200  // TODO FIXME somehow we have to transfer the striding information of the permuted domain/range maps.
201  // That is probably something for an external permutation factory
202  // if (originalP->IsView("stridedMaps"))
203  // rebalancedP->CreateView("stridedMaps", originalP);
205  if(!rebalancedP.is_null()) {std::ostringstream oss; oss << "P_" << coarseLevel.GetLevelID(); rebalancedP->setObjectLabel(oss.str());}
206  Set(coarseLevel, "P", rebalancedP);
207 
208  if (IsPrint(Statistics2))
209  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*rebalancedP, "P (rebalanced)", params);
210  }
211  }
212 
213  if (importer.is_null()) {
214  if (IsAvailable(coarseLevel, "Nullspace"))
215  Set(coarseLevel, "Nullspace", Get<RCP<MultiVector> >(coarseLevel, "Nullspace"));
216 
217  if (pL.isParameter("Coordinates") && pL.get< RCP<const FactoryBase> >("Coordinates") != Teuchos::null)
218  if (IsAvailable(coarseLevel, "Coordinates"))
219  Set(coarseLevel, "Coordinates", Get< RCP<xdMV> >(coarseLevel, "Coordinates"));
220 
221  return;
222  }
223 
224  if (pL.isParameter("Coordinates") &&
225  pL.get< RCP<const FactoryBase> >("Coordinates") != Teuchos::null &&
226  IsAvailable(coarseLevel, "Coordinates")) {
227  RCP<xdMV> coords = Get<RCP<xdMV> >(coarseLevel, "Coordinates");
228 
229  // This line must be after the Get call
230  SubFactoryMonitor subM(*this, "Rebalancing coordinates", coarseLevel);
231 
232  LO nodeNumElts = coords->getMap()->getNodeNumElements();
233 
234  // If a process has no matrix rows, then we can't calculate blocksize using the formula below.
235  LO myBlkSize = 0, blkSize = 0;
236  if (nodeNumElts > 0)
237  myBlkSize = importer->getSourceMap()->getNodeNumElements() / nodeNumElts;
238  MueLu_maxAll(coords->getMap()->getComm(), myBlkSize, blkSize);
239 
240  RCP<const Import> coordImporter;
241  if (blkSize == 1) {
242  coordImporter = importer;
243 
244  } else {
245  // NOTE: there is an implicit assumption here: we assume that dof any node are enumerated consequently
246  // Proper fix would require using decomposition similar to how we construct importer in the
247  // RepartitionFactory
248  RCP<const Map> origMap = coords->getMap();
249  GO indexBase = origMap->getIndexBase();
250 
251  ArrayView<const GO> OEntries = importer->getTargetMap()->getNodeElementList();
252  LO numEntries = OEntries.size()/blkSize;
253  ArrayRCP<GO> Entries(numEntries);
254  for (LO i = 0; i < numEntries; i++)
255  Entries[i] = (OEntries[i*blkSize]-indexBase)/blkSize + indexBase;
256 
257  RCP<const Map> targetMap = MapFactory::Build(origMap->lib(), origMap->getGlobalNumElements(), Entries(), indexBase, origMap->getComm());
258  coordImporter = ImportFactory::Build(origMap, targetMap);
259  }
260 
261  RCP<xdMV> permutedCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO>::Build(coordImporter->getTargetMap(), coords->getNumVectors());
262  permutedCoords->doImport(*coords, *coordImporter, Xpetra::INSERT);
263 
264  if (pL.isParameter("repartition: use subcommunicators") == true && pL.get<bool>("repartition: use subcommunicators") == true)
265  permutedCoords->replaceMap(permutedCoords->getMap()->removeEmptyProcesses());
266 
267  if (permutedCoords->getMap() == Teuchos::null)
268  permutedCoords = Teuchos::null;
269 
270  Set(coarseLevel, "Coordinates", permutedCoords);
271 
272  std::string fileName = "rebalanced_coordinates_level_" + toString(coarseLevel.GetLevelID()) + ".m";
273  if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd && permutedCoords->getMap() != Teuchos::null)
274  Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO>::Write(fileName, *permutedCoords);
275  }
276 
277  if (IsAvailable(coarseLevel, "Nullspace")) {
278  RCP<MultiVector> nullspace = Get< RCP<MultiVector> >(coarseLevel, "Nullspace");
279 
280  // This line must be after the Get call
281  SubFactoryMonitor subM(*this, "Rebalancing nullspace", coarseLevel);
282 
283  RCP<MultiVector> permutedNullspace = MultiVectorFactory::Build(importer->getTargetMap(), nullspace->getNumVectors());
284  permutedNullspace->doImport(*nullspace, *importer, Xpetra::INSERT);
285 
286  if (pL.get<bool>("repartition: use subcommunicators") == true)
287  permutedNullspace->replaceMap(permutedNullspace->getMap()->removeEmptyProcesses());
288 
289  if (permutedNullspace->getMap() == Teuchos::null)
290  permutedNullspace = Teuchos::null;
291 
292  Set(coarseLevel, "Nullspace", permutedNullspace);
293  }
294 
295  } else {
296  if (pL.get<bool>("transpose: use implicit") == false) {
297  RCP<Matrix> originalR = Get< RCP<Matrix> >(coarseLevel, "R");
298 
299  SubFactoryMonitor m2(*this, "Rebalancing restrictor", coarseLevel);
300 
301  if (implicit || importer.is_null()) {
302  GetOStream(Runtime0) << "Using original restrictor" << std::endl;
303  Set(coarseLevel, "R", originalR);
304 
305  } else {
306  RCP<Matrix> rebalancedR;
307  {
308  SubFactoryMonitor subM(*this, "Rebalancing restriction -- fusedImport", coarseLevel);
309 
310  RCP<Map> dummy; // meaning: use originalR's domain map.
311  Teuchos::ParameterList listLabel;
312  listLabel.set("Timer Label","MueLu::RebalanceR-" + Teuchos::toString(coarseLevel.GetLevelID()));
313  rebalancedR = MatrixFactory::Build(originalR, *importer, dummy, importer->getTargetMap(),Teuchos::rcp(&listLabel,false));
314  }
315  if(!rebalancedR.is_null()) {std::ostringstream oss; oss << "R_" << coarseLevel.GetLevelID(); rebalancedR->setObjectLabel(oss.str());}
316  Set(coarseLevel, "R", rebalancedR);
317 
319  // TODO FIXME somehow we have to transfer the striding information of the permuted domain/range maps.
320  // That is probably something for an external permutation factory
321  // if (originalR->IsView("stridedMaps"))
322  // rebalancedR->CreateView("stridedMaps", originalR);
324 
325  if (IsPrint(Statistics2))
326  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*rebalancedR, "R (rebalanced)", params);
327  }
328  }
329  }
330  }
331 
332 } // namespace MueLu
333 
334 #endif // MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
Exception indicating invalid cast attempted.
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-&gt;Get&lt; RCP&lt;Matrix&gt; &gt;(&quot;A&quot;, factory) if factory == NULL =&gt; use default factory.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
#define MueLu_maxAll(rcpComm, in, out)
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.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
size_type size() const
One-liner description of what is happening.
static const NoFactory * get()
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Print even more statistics.
bool isParameter(const std::string &name) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
#define SET_VALID_ENTRY(name)
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:76
std::string toString(const T &t)
bool is_null() const