MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_CloneRepartitionInterface_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 PACKAGES_MUELU_SRC_REBALANCING_MUELU_CLONEREPARTITIONINTERFACE_DEF_HPP_
11 #define PACKAGES_MUELU_SRC_REBALANCING_MUELU_CLONEREPARTITIONINTERFACE_DEF_HPP_
12 
13 #include <Xpetra_MultiVectorFactory.hpp>
14 #include <Xpetra_VectorFactory.hpp>
15 
17 
18 #include "MueLu_Level.hpp"
19 #include "MueLu_Exceptions.hpp"
20 #include "MueLu_Monitor.hpp"
21 
22 namespace MueLu {
23 
24 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
26  Teuchos::RCP<ParameterList> validParamList = rcp(new ParameterList());
27  validParamList->set<Teuchos::RCP<const FactoryBase> >("A", Teuchos::null, "Factory of the matrix A");
28  validParamList->set<Teuchos::RCP<const FactoryBase> >("number of partitions", Teuchos::null, "Instance of RepartitionHeuristicFactory.");
29  validParamList->set<Teuchos::RCP<const FactoryBase> >("Partition", Teuchos::null, "Factory generating the Partition array (e.g. ZoltanInterface)");
30 
31  return validParamList;
32 }
33 
34 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
36  Input(currentLevel, "A");
37  Input(currentLevel, "Partition");
38 } // DeclareInput()
39 
40 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
42  FactoryMonitor m(*this, "Build", currentLevel);
43  currentLevel.print(GetOStream(Statistics0, 0));
44  // extract blocked operator A from current level
45  Teuchos::RCP<Matrix> A = Get<Teuchos::RCP<Matrix> >(currentLevel, "A");
46  Teuchos::RCP<const Teuchos::Comm<int> > comm = A->getRowMap()->getComm();
47 
48  // number of Partitions only used for a shortcut.
49  if (currentLevel.IsAvailable("number of partitions")) {
50  GetOStream(Warnings0) << "Using user-provided \"number of partitions\", the performance is unknown" << std::endl;
51  }
52 
53  // ======================================================================================================
54  // Construct decomposition vector
55  // ======================================================================================================
56  RCP<GOVector> decomposition = Teuchos::null;
57 
58  // extract decomposition vector
59  decomposition = Get<RCP<GOVector> >(currentLevel, "Partition");
60  ArrayRCP<const GO> decompEntries = decomposition->getData(0);
61 
62  if (decomposition.is_null()) {
63  GetOStream(Warnings0) << "No repartitioning necessary: partitions were left unchanged by the repartitioner" << std::endl;
64  Set<RCP<const Import> >(currentLevel, "Importer", Teuchos::null);
65  return;
66  }
67 
68  // create new decomposition vector
70  ArrayRCP<GO> retDecompEntries = ret->getDataNonConst(0);
71 
72  // block size of output vector
73  LocalOrdinal blkSize = 1;
74 
75  // check for blocking/striding information
76  if (A->IsView("stridedMaps") &&
77  Teuchos::rcp_dynamic_cast<const StridedMap>(A->getRowMap("stridedMaps")) != Teuchos::null) {
78  Xpetra::viewLabel_t oldView = A->SwitchToView("stridedMaps"); // note: "stridedMaps are always non-overlapping (correspond to range and domain maps!)
79  RCP<const StridedMap> strMap = Teuchos::rcp_dynamic_cast<const StridedMap>(A->getRowMap());
80  TEUCHOS_TEST_FOR_EXCEPTION(strMap == Teuchos::null, Exceptions::BadCast, "MueLu::CloneRepartitionInterface::Build: cast to strided row map failed.");
81  LocalOrdinal stridedBlock = strMap->getStridedBlockId();
82  if (stridedBlock == -1)
83  blkSize = strMap->getFixedBlockSize();
84  else {
85  std::vector<size_t> strInfo = strMap->getStridingData();
86  blkSize = strInfo[stridedBlock];
87  }
88  oldView = A->SwitchToView(oldView);
89  GetOStream(Statistics1) << "CloneRepartitionInterface::Build():"
90  << " found blockdim=" << blkSize << " from strided maps." << std::endl;
91  } else {
92  GetOStream(Statistics1) << "CloneRepartitionInterface::Build(): no striding information available. Use blockdim=" << blkSize << " (DofsPerNode)." << std::endl;
93  blkSize = A->GetFixedBlockSize();
94  }
95 
96  // plausibility check!
97  size_t inLocalLength = decomposition->getLocalLength();
98  size_t outLocalLength = A->getRowMap()->getLocalNumElements();
99 
100  // only for non-strided maps
101  size_t numLocalNodes = outLocalLength / blkSize;
102  TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::as<size_t>(outLocalLength % blkSize) != 0, MueLu::Exceptions::RuntimeError, "CloneRepartitionInterface: inconsistent number of local DOFs (" << outLocalLength << ") and degrees of freedoms (" << blkSize << ")");
103 
104  if (numLocalNodes > 0) {
105  TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::as<size_t>(inLocalLength % numLocalNodes) != 0, MueLu::Exceptions::RuntimeError, "CloneRepartitionInterface: inconsistent number of local DOFs (" << inLocalLength << ") and number of local nodes (" << numLocalNodes << ")");
106  LocalOrdinal inBlkSize = Teuchos::as<LocalOrdinal>(inLocalLength / numLocalNodes);
107  // TEUCHOS_TEST_FOR_EXCEPTION(blkSize != inBlkSize, MueLu::Exceptions::RuntimeError,"CloneRepartitionInterface: input block size = " << inBlkSize << " outpub block size = " << blkSize << ". They should be the same.");
108 
109  for (LO i = 0; i < Teuchos::as<LO>(numLocalNodes); i++) {
110  for (LO j = 0; j < blkSize; j++) {
111  retDecompEntries[i * blkSize + j] = Teuchos::as<GO>(decompEntries[i * inBlkSize]);
112  }
113  }
114  } // end if numLocalNodes > 0
115  Set(currentLevel, "Partition", ret);
116 } // Build()
117 
118 } // namespace MueLu
119 
120 #endif /* PACKAGES_MUELU_SRC_REBALANCING_MUELU_CLONEREPARTITIONINTERFACE_DEF_HPP_ */
Important warning messages (one line)
Exception indicating invalid cast attempted.
MueLu::DefaultLocalOrdinal LocalOrdinal
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Print more statistics.
LocalOrdinal LO
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
void DeclareInput(Level &level) const
Specifies the data that this class needs, and the factories that generate that data.
Print statistics that do not involve significant additional computation.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
std::string viewLabel_t
static RCP< Vector > Build(const Teuchos::RCP< const Map > &map, bool zeroOut=true)
void print(std::ostream &out, const VerbLevel verbLevel=Default) const
Printing method.
Exception throws to report errors in the internal logical of the program.
void Build(Level &level) const
Build an object with this factory.
bool is_null() const