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 // ***********************************************************************
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_COORDINATESTRANSFER_FACTORY_DEF_HPP
47 #define MUELU_COORDINATESTRANSFER_FACTORY_DEF_HPP
48 
49 #include "Xpetra_ImportFactory.hpp"
50 #include "Xpetra_MultiVectorFactory.hpp"
51 #include "Xpetra_MapFactory.hpp"
52 #include "Xpetra_IO.hpp"
53 
54 #include "MueLu_Aggregates.hpp"
56 #include "MueLu_Utilities.hpp"
57 
58 #include "MueLu_Level.hpp"
59 #include "MueLu_Monitor.hpp"
60 
61 namespace MueLu {
62 
63 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
65  RCP<ParameterList> validParamList = rcp(new ParameterList());
66 
67  validParamList->set<RCP<const FactoryBase> >("Coordinates", Teuchos::null, "Factory for coordinates generation");
68  validParamList->set<RCP<const FactoryBase> >("Aggregates", Teuchos::null, "Factory for coordinates generation");
69  validParamList->set<RCP<const FactoryBase> >("CoarseMap", Teuchos::null, "Generating factory of the coarse map");
70  validParamList->set<bool>("structured aggregation", false, "Flag specifying that the geometric data is transferred for StructuredAggregationFactory");
71  validParamList->set<bool>("aggregation coupled", false, "Flag specifying if the aggregation algorithm was used in coupled mode.");
72  validParamList->set<bool>("Geometric", false, "Flag specifying that the coordinates are transferred for GeneralGeometricPFactory");
73  validParamList->set<RCP<const FactoryBase> >("coarseCoordinates", Teuchos::null, "Factory for coarse coordinates generation");
74  validParamList->set<RCP<const FactoryBase> >("gCoarseNodesPerDim", Teuchos::null, "Factory providing the global number of nodes per spatial dimensions of the mesh");
75  validParamList->set<RCP<const FactoryBase> >("lCoarseNodesPerDim", Teuchos::null, "Factory providing the local number of nodes per spatial dimensions of the mesh");
76  validParamList->set<RCP<const FactoryBase> >("numDimensions", Teuchos::null, "Factory providing the number of spatial dimensions of the mesh");
77  validParamList->set<int>("write start", -1, "first level at which coordinates should be written to file");
78  validParamList->set<int>("write end", -1, "last level at which coordinates should be written to file");
79  validParamList->set<bool>("hybrid aggregation", false, "Flag specifying that hybrid aggregation data is transfered for HybridAggregationFactory");
80  validParamList->set<RCP<const FactoryBase> >("aggregationRegionTypeCoarse", Teuchos::null, "Factory indicating what aggregation type is to be used on the coarse level of the region");
81  validParamList->set<bool>("interface aggregation", false, "Flag specifying that interface aggregation data is transfered for HybridAggregationFactory");
82  validParamList->set<RCP<const FactoryBase> >("coarseInterfacesDimensions", Teuchos::null, "Factory providing coarseInterfacesDimensions");
83  validParamList->set<RCP<const FactoryBase> >("nodeOnCoarseInterface", Teuchos::null, "Factory providing nodeOnCoarseInterface");
84 
85  return validParamList;
86 }
87 
88 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
90  static bool isAvailableCoords = false;
91 
92  const ParameterList& pL = GetParameterList();
93  if (pL.get<bool>("structured aggregation") == true) {
94  if (pL.get<bool>("aggregation coupled") == true) {
95  Input(fineLevel, "gCoarseNodesPerDim");
96  }
97  Input(fineLevel, "lCoarseNodesPerDim");
98  Input(fineLevel, "numDimensions");
99  } else if (pL.get<bool>("Geometric") == true) {
100  Input(coarseLevel, "coarseCoordinates");
101  Input(coarseLevel, "gCoarseNodesPerDim");
102  Input(coarseLevel, "lCoarseNodesPerDim");
103  } else if (pL.get<bool>("hybrid aggregation") == true) {
104  Input(fineLevel, "aggregationRegionTypeCoarse");
105  Input(fineLevel, "lCoarseNodesPerDim");
106  Input(fineLevel, "numDimensions");
107  if (pL.get<bool>("interface aggregation") == true) {
108  Input(fineLevel, "coarseInterfacesDimensions");
109  Input(fineLevel, "nodeOnCoarseInterface");
110  }
111  } else {
112  if (coarseLevel.GetRequestMode() == Level::REQUEST)
113  isAvailableCoords = coarseLevel.IsAvailable("Coordinates", this);
114 
115  if (isAvailableCoords == false) {
116  Input(fineLevel, "Coordinates");
117  Input(fineLevel, "Aggregates");
118  Input(fineLevel, "CoarseMap");
119  }
120  }
121 }
122 
123 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
125  FactoryMonitor m(*this, "Build", coarseLevel);
126 
128 
129  GetOStream(Runtime0) << "Transferring coordinates" << std::endl;
130 
131  int numDimensions;
132  RCP<xdMV> coarseCoords;
133  RCP<xdMV> fineCoords;
134  Array<GO> gCoarseNodesPerDir;
135  Array<LO> lCoarseNodesPerDir;
136 
137  const ParameterList& pL = GetParameterList();
138 
139  if (pL.get<bool>("hybrid aggregation") == true) {
140  std::string regionType = Get<std::string>(fineLevel, "aggregationRegionTypeCoarse");
141  numDimensions = Get<int>(fineLevel, "numDimensions");
142  lCoarseNodesPerDir = Get<Array<LO> >(fineLevel, "lCoarseNodesPerDim");
143  Set<std::string>(coarseLevel, "aggregationRegionType", regionType);
144  Set<int>(coarseLevel, "numDimensions", numDimensions);
145  Set<Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
146 
147  if ((pL.get<bool>("interface aggregation") == true) && (regionType == "uncoupled")) {
148  Array<LO> coarseInterfacesDimensions = Get<Array<LO> >(fineLevel, "coarseInterfacesDimensions");
149  Array<LO> nodeOnCoarseInterface = Get<Array<LO> >(fineLevel, "nodeOnCoarseInterface");
150  Set<Array<LO> >(coarseLevel, "interfacesDimensions", coarseInterfacesDimensions);
151  Set<Array<LO> >(coarseLevel, "nodeOnInterface", nodeOnCoarseInterface);
152  }
153 
154  } else if (pL.get<bool>("structured aggregation") == true) {
155  if (pL.get<bool>("aggregation coupled") == true) {
156  gCoarseNodesPerDir = Get<Array<GO> >(fineLevel, "gCoarseNodesPerDim");
157  Set<Array<GO> >(coarseLevel, "gNodesPerDim", gCoarseNodesPerDir);
158  }
159  lCoarseNodesPerDir = Get<Array<LO> >(fineLevel, "lCoarseNodesPerDim");
160  Set<Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
161  numDimensions = Get<int>(fineLevel, "numDimensions");
162  Set<int>(coarseLevel, "numDimensions", numDimensions);
163 
164  } else if (pL.get<bool>("Geometric") == true) {
165  coarseCoords = Get<RCP<xdMV> >(coarseLevel, "coarseCoordinates");
166  gCoarseNodesPerDir = Get<Array<GO> >(coarseLevel, "gCoarseNodesPerDim");
167  lCoarseNodesPerDir = Get<Array<LO> >(coarseLevel, "lCoarseNodesPerDim");
168  Set<Array<GO> >(coarseLevel, "gNodesPerDim", gCoarseNodesPerDir);
169  Set<Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
170 
171  Set<RCP<xdMV> >(coarseLevel, "Coordinates", coarseCoords);
172 
173  } else {
174  if (coarseLevel.IsAvailable("Coordinates", this)) {
175  GetOStream(Runtime0) << "Reusing coordinates" << std::endl;
176  return;
177  }
178 
179  fineCoords = Get<RCP<xdMV> >(fineLevel, "Coordinates");
180  RCP<const Map> coarseMap = Get<RCP<const Map> >(fineLevel, "CoarseMap");
181 
182  // coarseMap is being used to set up the domain map of tentative P, and therefore, the row map of Ac
183  // Therefore, if we amalgamate coarseMap, logical nodes in the coordinates vector would correspond to
184  // logical blocks in the matrix
185  LO blkSize = 1;
186  if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null)
187  blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap)->getFixedBlockSize();
188 
189  RCP<const Map> coarseCoordMap;
190  RCP<const Map> uniqueMap = fineCoords->getMap();
191  if (blkSize > 1) {
192  // If the block size is greater than one, we need to create a coarse coordinate map
193  // FIXME: The amalgamation should really be done on device.
194  GO indexBase = coarseMap->getIndexBase();
195  ArrayView<const GO> elementAList = coarseMap->getLocalElementList();
196  size_t numElements = elementAList.size() / blkSize;
197  Array<GO> elementList(numElements);
198 
199  // Amalgamate the map
200  for (LO i = 0; i < Teuchos::as<LO>(numElements); i++)
201  elementList[i] = (elementAList[i * blkSize] - indexBase) / blkSize + indexBase;
202 
203  {
204  SubFactoryMonitor sfm(*this, "MapFactory: coarseCoordMap", fineLevel);
205  coarseCoordMap = MapFactory ::Build(coarseMap->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), elementList, indexBase, coarseMap->getComm());
206  }
207  } else {
208  // If the block size is one, we can just use the coarse map for coordinates
209  coarseCoordMap = coarseMap;
210  }
211 
212  // Build the coarseCoords MultiVector
213  coarseCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Build(coarseCoordMap, fineCoords->getNumVectors());
214 
215  RCP<Aggregates> aggregates;
216  bool aggregatesCrossProcessors;
217  aggregates = Get<RCP<Aggregates> >(fineLevel, "Aggregates");
218  aggregatesCrossProcessors = aggregates->AggregatesCrossProcessors();
219 
220  // Create overlapped fine coordinates to reduce global communication
221  RCP<xdMV> ghostedCoords = fineCoords;
222  if (aggregatesCrossProcessors) {
223  RCP<const Map> nonUniqueMap = aggregates->GetMap();
224  RCP<const Import> importer = ImportFactory::Build(uniqueMap, nonUniqueMap);
225 
226  ghostedCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Build(nonUniqueMap, fineCoords->getNumVectors());
227  ghostedCoords->doImport(*fineCoords, *importer, Xpetra::INSERT);
228  }
229 
230  // The good news is that this graph has already been constructed for the
231  // TentativePFactory and was cached in Aggregates. So this is a no-op.
232  auto aggGraph = aggregates->GetGraph();
233  auto numAggs = aggGraph.numRows();
234 
235  auto fineCoordsView = ghostedCoords->getDeviceLocalView(Xpetra::Access::ReadOnly);
236  auto coarseCoordsView = coarseCoords->getDeviceLocalView(Xpetra::Access::OverwriteAll);
237 
238  // Fill in coarse coordinates
239  {
240  SubFactoryMonitor m2(*this, "AverageCoords", coarseLevel);
241 
242  const auto dim = ghostedCoords->getNumVectors();
243 
244  typename AppendTrait<decltype(fineCoordsView), Kokkos::RandomAccess>::type fineCoordsRandomView = fineCoordsView;
245  for (size_t j = 0; j < dim; j++) {
246  Kokkos::parallel_for(
247  "MueLu:CoordinatesTransferF:Build:coord", Kokkos::RangePolicy<local_ordinal_type, execution_space>(0, numAggs),
248  KOKKOS_LAMBDA(const LO i) {
249  // A row in this graph represents all node ids in the aggregate
250  // Therefore, averaging is very easy
251 
252  auto aggregate = aggGraph.rowConst(i);
253 
254  typename Teuchos::ScalarTraits<Scalar>::magnitudeType sum = 0.0; // do not use Scalar here (Stokhos)
255  for (size_t colID = 0; colID < static_cast<size_t>(aggregate.length); colID++)
256  sum += fineCoordsRandomView(aggregate(colID), j);
257 
258  coarseCoordsView(i, j) = sum / aggregate.length;
259  });
260  }
261  }
262 
263  Set<RCP<xdMV> >(coarseLevel, "Coordinates", coarseCoords);
264  }
265 
266  int writeStart = pL.get<int>("write start"), writeEnd = pL.get<int>("write end");
267  if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd) {
268  std::ostringstream buf;
269  buf << fineLevel.GetLevelID();
270  std::string fileName = "coordinates_before_rebalance_level_" + buf.str() + ".m";
271  Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Write(fileName, *fineCoords);
272  }
273  if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd) {
274  std::ostringstream buf;
275  buf << coarseLevel.GetLevelID();
276  std::string fileName = "coordinates_before_rebalance_level_" + buf.str() + ".m";
277  Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Write(fileName, *coarseCoords);
278  }
279 }
280 
281 } // namespace MueLu
282 
283 #endif // MUELU_COORDINATESTRANSFER_FACTORY_DEF_HPP
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)
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.
size_type size() const
LocalOrdinal LO
One-liner description of what is happening.
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:99
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
Node NO
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:76
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.