Loading [MathJax]/extensions/tex2jax.js
MueLu  Version of the Day
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_CoordinatesTransferFactory_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_COORDINATESTRANSFER_FACTORY_DEF_HPP
11 #define MUELU_COORDINATESTRANSFER_FACTORY_DEF_HPP
12 
13 #include "Xpetra_ImportFactory.hpp"
14 #include "Xpetra_MultiVectorFactory.hpp"
15 #include "Xpetra_MapFactory.hpp"
16 #include "Xpetra_IO.hpp"
17 
18 #include "MueLu_Aggregates.hpp"
19 #include "MueLu_AmalgamationFactory.hpp"
21 #include "MueLu_Utilities.hpp"
22 
23 #include "MueLu_Level.hpp"
24 #include "MueLu_Monitor.hpp"
25 
26 namespace MueLu {
27 
28 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
30 
31 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
33 
34 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
36  RCP<ParameterList> validParamList = rcp(new ParameterList());
37 
38  validParamList->set<RCP<const FactoryBase> >("Coordinates", Teuchos::null, "Factory for coordinates generation");
39  validParamList->set<RCP<const FactoryBase> >("Aggregates", Teuchos::null, "Factory for coordinates generation");
40  validParamList->set<RCP<const FactoryBase> >("CoarseMap", Teuchos::null, "Generating factory of the coarse map");
41  validParamList->set<bool>("structured aggregation", false, "Flag specifying that the geometric data is transferred for StructuredAggregationFactory");
42  validParamList->set<bool>("aggregation coupled", false, "Flag specifying if the aggregation algorithm was used in coupled mode.");
43  validParamList->set<bool>("Geometric", false, "Flag specifying that the coordinates are transferred for GeneralGeometricPFactory");
44  validParamList->set<RCP<const FactoryBase> >("coarseCoordinates", Teuchos::null, "Factory for coarse coordinates generation");
45  validParamList->set<RCP<const FactoryBase> >("gCoarseNodesPerDim", Teuchos::null, "Factory providing the global number of nodes per spatial dimensions of the mesh");
46  validParamList->set<RCP<const FactoryBase> >("lCoarseNodesPerDim", Teuchos::null, "Factory providing the local number of nodes per spatial dimensions of the mesh");
47  validParamList->set<RCP<const FactoryBase> >("numDimensions", Teuchos::null, "Factory providing the number of spatial dimensions of the mesh");
48  validParamList->set<int>("write start", -1, "first level at which coordinates should be written to file");
49  validParamList->set<int>("write end", -1, "last level at which coordinates should be written to file");
50  validParamList->set<bool>("hybrid aggregation", false, "Flag specifying that hybrid aggregation data is transfered for HybridAggregationFactory");
51  validParamList->set<RCP<const FactoryBase> >("aggregationRegionTypeCoarse", Teuchos::null, "Factory indicating what aggregation type is to be used on the coarse level of the region");
52  validParamList->set<bool>("interface aggregation", false, "Flag specifying that interface aggregation data is transfered for HybridAggregationFactory");
53  validParamList->set<RCP<const FactoryBase> >("coarseInterfacesDimensions", Teuchos::null, "Factory providing coarseInterfacesDimensions");
54  validParamList->set<RCP<const FactoryBase> >("nodeOnCoarseInterface", Teuchos::null, "Factory providing nodeOnCoarseInterface");
55 
56  return validParamList;
57 }
58 
59 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
61  static bool isAvailableCoords = false;
62 
63  const ParameterList& pL = GetParameterList();
64  if (pL.get<bool>("structured aggregation") == true) {
65  if (pL.get<bool>("aggregation coupled") == true) {
66  Input(fineLevel, "gCoarseNodesPerDim");
67  }
68  Input(fineLevel, "lCoarseNodesPerDim");
69  Input(fineLevel, "numDimensions");
70  } else if (pL.get<bool>("Geometric") == true) {
71  Input(coarseLevel, "coarseCoordinates");
72  Input(coarseLevel, "gCoarseNodesPerDim");
73  Input(coarseLevel, "lCoarseNodesPerDim");
74  } else if (pL.get<bool>("hybrid aggregation") == true) {
75  Input(fineLevel, "aggregationRegionTypeCoarse");
76  Input(fineLevel, "lCoarseNodesPerDim");
77  Input(fineLevel, "numDimensions");
78  if (pL.get<bool>("interface aggregation") == true) {
79  Input(fineLevel, "coarseInterfacesDimensions");
80  Input(fineLevel, "nodeOnCoarseInterface");
81  }
82  } else {
83  if (coarseLevel.GetRequestMode() == Level::REQUEST)
84  isAvailableCoords = coarseLevel.IsAvailable("Coordinates", this);
85 
86  if (isAvailableCoords == false) {
87  Input(fineLevel, "Coordinates");
88  Input(fineLevel, "Aggregates");
89  Input(fineLevel, "CoarseMap");
90  }
91  }
92 }
93 
94 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
96  FactoryMonitor m(*this, "Build", coarseLevel);
97 
99 
100  GetOStream(Runtime0) << "Transferring coordinates" << std::endl;
101 
102  int numDimensions;
103  RCP<xdMV> coarseCoords;
104  RCP<xdMV> fineCoords;
105  Array<GO> gCoarseNodesPerDir;
106  Array<LO> lCoarseNodesPerDir;
107 
108  const ParameterList& pL = GetParameterList();
109 
110  if (pL.get<bool>("hybrid aggregation") == true) {
111  std::string regionType = Get<std::string>(fineLevel, "aggregationRegionTypeCoarse");
112  numDimensions = Get<int>(fineLevel, "numDimensions");
113  lCoarseNodesPerDir = Get<Array<LO> >(fineLevel, "lCoarseNodesPerDim");
114  Set<std::string>(coarseLevel, "aggregationRegionType", regionType);
115  Set<int>(coarseLevel, "numDimensions", numDimensions);
116  Set<Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
117 
118  if ((pL.get<bool>("interface aggregation") == true) && (regionType == "uncoupled")) {
119  Array<LO> coarseInterfacesDimensions = Get<Array<LO> >(fineLevel, "coarseInterfacesDimensions");
120  Array<LO> nodeOnCoarseInterface = Get<Array<LO> >(fineLevel, "nodeOnCoarseInterface");
121  Set<Array<LO> >(coarseLevel, "interfacesDimensions", coarseInterfacesDimensions);
122  Set<Array<LO> >(coarseLevel, "nodeOnInterface", nodeOnCoarseInterface);
123  }
124 
125  } else if (pL.get<bool>("structured aggregation") == true) {
126  if (pL.get<bool>("aggregation coupled") == true) {
127  gCoarseNodesPerDir = Get<Array<GO> >(fineLevel, "gCoarseNodesPerDim");
128  Set<Array<GO> >(coarseLevel, "gNodesPerDim", gCoarseNodesPerDir);
129  }
130  lCoarseNodesPerDir = Get<Array<LO> >(fineLevel, "lCoarseNodesPerDim");
131  Set<Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
132  numDimensions = Get<int>(fineLevel, "numDimensions");
133  Set<int>(coarseLevel, "numDimensions", numDimensions);
134 
135  } else if (pL.get<bool>("Geometric") == true) {
136  coarseCoords = Get<RCP<xdMV> >(coarseLevel, "coarseCoordinates");
137  gCoarseNodesPerDir = Get<Array<GO> >(coarseLevel, "gCoarseNodesPerDim");
138  lCoarseNodesPerDir = Get<Array<LO> >(coarseLevel, "lCoarseNodesPerDim");
139  Set<Array<GO> >(coarseLevel, "gNodesPerDim", gCoarseNodesPerDir);
140  Set<Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
141 
142  Set<RCP<xdMV> >(coarseLevel, "Coordinates", coarseCoords);
143 
144  } else {
145  if (coarseLevel.IsAvailable("Coordinates", this)) {
146  GetOStream(Runtime0) << "Reusing coordinates" << std::endl;
147  return;
148  }
149 
150  fineCoords = Get<RCP<xdMV> >(fineLevel, "Coordinates");
151  RCP<const Map> coarseMap = Get<RCP<const Map> >(fineLevel, "CoarseMap");
152 
153  RCP<const Map> coarseCoordMap;
154  RCP<const Map> uniqueMap = fineCoords->getMap();
155 
156  // coarseMap is being used to set up the domain map of tentative P, and therefore, the row map of Ac
157  // Therefore, if we amalgamate coarseMap, logical nodes in the coordinates vector would correspond to
158  // logical blocks in the matrix
159  LO blkSize = 1;
160  if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null)
161  blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap)->getFixedBlockSize();
162 
163  if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null)
164  blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap)->getFixedBlockSize();
165 
166  if (blkSize == 1) {
167  // Scalar system
168  // No amalgamation required, we can use the coarseMap
169  coarseCoordMap = coarseMap;
170  } else {
171  // Vector system
172  AmalgamationFactory::AmalgamateMap(rcp_dynamic_cast<const StridedMap>(coarseMap), coarseCoordMap);
173  }
174 
175  // Build the coarseCoords MultiVector
176  coarseCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Build(coarseCoordMap, fineCoords->getNumVectors());
177 
178  RCP<Aggregates> aggregates;
179  bool aggregatesCrossProcessors;
180  aggregates = Get<RCP<Aggregates> >(fineLevel, "Aggregates");
181  aggregatesCrossProcessors = aggregates->AggregatesCrossProcessors();
182 
183  // Create overlapped fine coordinates to reduce global communication
184  RCP<xdMV> ghostedCoords = fineCoords;
185  if (aggregatesCrossProcessors) {
186  RCP<const Map> nonUniqueMap = aggregates->GetMap();
187  RCP<const Import> importer = ImportFactory::Build(uniqueMap, nonUniqueMap);
188 
189  ghostedCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Build(nonUniqueMap, fineCoords->getNumVectors());
190  ghostedCoords->doImport(*fineCoords, *importer, Xpetra::INSERT);
191  }
192 
193  // The good news is that this graph has already been constructed for the
194  // TentativePFactory and was cached in Aggregates. So this is a no-op.
195  auto aggGraph = aggregates->GetGraph();
196  auto numAggs = aggGraph.numRows();
197 
198  auto fineCoordsView = ghostedCoords->getDeviceLocalView(Xpetra::Access::ReadOnly);
199  auto coarseCoordsView = coarseCoords->getDeviceLocalView(Xpetra::Access::OverwriteAll);
200 
201  // Fill in coarse coordinates
202  {
203  SubFactoryMonitor m2(*this, "AverageCoords", coarseLevel);
204 
205  const auto dim = ghostedCoords->getNumVectors();
206 
207  typename AppendTrait<decltype(fineCoordsView), Kokkos::RandomAccess>::type fineCoordsRandomView = fineCoordsView;
208  for (size_t j = 0; j < dim; j++) {
209  Kokkos::parallel_for(
210  "MueLu:CoordinatesTransferF:Build:coord", Kokkos::RangePolicy<local_ordinal_type, execution_space>(0, numAggs),
211  KOKKOS_LAMBDA(const LO i) {
212  // A row in this graph represents all node ids in the aggregate
213  // Therefore, averaging is very easy
214 
215  auto aggregate = aggGraph.rowConst(i);
216 
217  typename Teuchos::ScalarTraits<Scalar>::magnitudeType sum = 0.0; // do not use Scalar here (Stokhos)
218  for (size_t colID = 0; colID < static_cast<size_t>(aggregate.length); colID++)
219  sum += fineCoordsRandomView(aggregate(colID), j);
220 
221  coarseCoordsView(i, j) = sum / aggregate.length;
222  });
223  }
224  }
225 
226  Set<RCP<xdMV> >(coarseLevel, "Coordinates", coarseCoords);
227  }
228 
229  int writeStart = pL.get<int>("write start"), writeEnd = pL.get<int>("write end");
230  if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd) {
231  std::ostringstream buf;
232  buf << fineLevel.GetLevelID();
233  std::string fileName = "coordinates_before_rebalance_level_" + buf.str() + ".m";
234  Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Write(fileName, *fineCoords);
235  }
236  if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd) {
237  std::ostringstream buf;
238  buf << coarseLevel.GetLevelID();
239  std::string fileName = "coordinates_before_rebalance_level_" + buf.str() + ".m";
240  Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Write(fileName, *coarseCoords);
241  }
242 }
243 
244 } // namespace MueLu
245 
246 #endif // MUELU_COORDINATESTRANSFER_FACTORY_DEF_HPP
CoordinatesTransferFactory()
Constructor.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
GlobalOrdinal GO
void DeclareInput(Level &finelevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
T & get(const std::string &name, T def_value)
Timer to be used in factories. Similar to Monitor but with additional timers.
LocalOrdinal LO
One-liner description of what is happening.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
virtual ~CoordinatesTransferFactory()
Destructor.
RequestMode GetRequestMode() const
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
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
static void AmalgamateMap(const Map &sourceMap, const Matrix &A, RCP< const Map > &amalgamatedMap, Array< LO > &translation)
Method to create merged map for systems of PDEs.
Node NO
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:51
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need&#39;s value has been saved.