MueLu  Version of the Day
MueLu_CoordinatesTransferFactory_kokkos_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_KOKKOS_DEF_HPP
47 #define MUELU_COORDINATESTRANSFER_FACTORY_KOKKOS_DEF_HPP
48 
50 
51 #include <Xpetra_ImportFactory.hpp>
52 #include <Xpetra_IO.hpp>
53 #include <Xpetra_MultiVectorFactory.hpp>
54 #include <Xpetra_MapFactory.hpp>
55 
56 #include "MueLu_Aggregates_kokkos.hpp"
57 #include "MueLu_CoarseMapFactory_kokkos.hpp"
58 #include "MueLu_Level.hpp"
59 #include "MueLu_Monitor.hpp"
60 #include "MueLu_Utilities_kokkos.hpp"
61 
62 namespace MueLu {
63 
64  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class DeviceType>
65  RCP<const ParameterList> CoordinatesTransferFactory_kokkos<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::GetValidParameterList() const {
66  RCP<ParameterList> validParamList = rcp(new ParameterList());
67 
68  validParamList->set< RCP<const FactoryBase> >("Coordinates", Teuchos::null, "Factory for coordinates generation");
69  validParamList->set< RCP<const FactoryBase> >("Aggregates", Teuchos::null, "Factory for coordinates generation");
70  validParamList->set< RCP<const FactoryBase> >("CoarseMap", Teuchos::null, "Generating factory of the coarse map");
71  validParamList->set< int > ("write start", -1, "First level at which coordinates should be written to file");
72  validParamList->set< int > ("write end", -1, "Last level at which coordinates should be written to file");
73 
74  return validParamList;
75  }
76 
77  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class DeviceType>
78  void CoordinatesTransferFactory_kokkos<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::DeclareInput(Level& fineLevel, Level& coarseLevel) const {
79  static bool isAvailableCoords = false;
80 
81  if (coarseLevel.GetRequestMode() == Level::REQUEST)
82  isAvailableCoords = coarseLevel.IsAvailable("Coordinates", this);
83 
84  if (isAvailableCoords == false) {
85  Input(fineLevel, "Coordinates");
86  Input(fineLevel, "Aggregates");
87  Input(fineLevel, "CoarseMap");
88  }
89  }
90 
91  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class DeviceType>
92  void CoordinatesTransferFactory_kokkos<Scalar, LocalOrdinal, GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::Build(Level& fineLevel, Level& coarseLevel) const {
93  FactoryMonitor m(*this, "Build", coarseLevel);
94 
95  typedef Xpetra::MultiVector<double,LO,GO,NO> doubleMultiVector;
96 
97  GetOStream(Runtime0) << "Transferring coordinates" << std::endl;
98 
99  if (coarseLevel.IsAvailable("Coordinates", this)) {
100  GetOStream(Runtime0) << "Reusing coordinates" << std::endl;
101  return;
102  }
103 
104  auto aggregates = Get<RCP<Aggregates_kokkos> >(fineLevel, "Aggregates");
105  auto fineCoords = Get<RCP<doubleMultiVector> >(fineLevel, "Coordinates");
106  auto coarseMap = Get<RCP<const Map> > (fineLevel, "CoarseMap");
107 
108  // coarseMap is being used to set up the domain map of tentative P, and therefore, the row map of Ac
109  // Therefore, if we amalgamate coarseMap, logical nodes in the coordinates vector would correspond to
110  // logical blocks in the matrix
111 
112  ArrayView<const GO> elementAList = coarseMap->getNodeElementList();
113  GO indexBase = coarseMap->getIndexBase();
114 
115  LO blkSize = 1;
116  if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null)
117  blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap)->getFixedBlockSize();
118 
119  Array<GO> elementList;
120  ArrayView<const GO> elementListView;
121  if (blkSize == 1) {
122  // Scalar system
123  // No amalgamation required
124  elementListView = elementAList;
125 
126  } else {
127  auto numElements = elementAList.size() / blkSize;
128 
129  elementList.resize(numElements);
130 
131  // Amalgamate the map
132  for (LO i = 0; i < Teuchos::as<LO>(numElements); i++)
133  elementList[i] = (elementAList[i*blkSize]-indexBase)/blkSize + indexBase;
134 
135  elementListView = elementList;
136  }
137 
138  auto uniqueMap = fineCoords->getMap();
139  auto coarseCoordMap = MapFactory::Build(coarseMap->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
140  elementListView, indexBase, coarseMap->getComm());
141  RCP<doubleMultiVector> coarseCoords = Xpetra::MultiVectorFactory<double,LO,GO,NO>::Build(coarseCoordMap, fineCoords->getNumVectors());
142 
143  // Create overlapped fine coordinates to reduce global communication
144  RCP<doubleMultiVector> ghostedCoords = fineCoords;
145  if (aggregates->AggregatesCrossProcessors()) {
146  auto nonUniqueMap = aggregates->GetMap();
147  auto importer = ImportFactory::Build(uniqueMap, nonUniqueMap);
148 
149  ghostedCoords = Xpetra::MultiVectorFactory<double,LO,GO,NO>::Build(nonUniqueMap, fineCoords->getNumVectors());
150  ghostedCoords->doImport(*fineCoords, *importer, Xpetra::INSERT);
151  }
152 
153  // The good new is that his graph has already been constructed for the
154  // TentativePFactory and was cached in Aggregates. So this is a no-op.
155  auto aggGraph = aggregates->GetGraph();
156  auto numAggs = aggGraph.numRows();
157 
158  auto fineCoordsView = fineCoords ->template getLocalView<DeviceType>();
159  auto coarseCoordsView = coarseCoords->template getLocalView<DeviceType>();
160 
161  // Fill in coarse coordinates
162  {
163  SubFactoryMonitor m2(*this, "AverageCoords", coarseLevel);
164 
165  const auto dim = fineCoords->getNumVectors();
166 
167  typename AppendTrait<decltype(fineCoordsView), Kokkos::RandomAccess>::type fineCoordsRandomView = fineCoordsView;
168  for (size_t j = 0; j < dim; j++) {
169  Kokkos::parallel_for("MueLu:CoordinatesTransferF:Build:coord", Kokkos::RangePolicy<local_ordinal_type, execution_space>(0, numAggs),
170  KOKKOS_LAMBDA(const LO i) {
171  // A row in this graph represents all node ids in the aggregate
172  // Therefore, averaging is very easy
173 
174  auto aggregate = aggGraph.rowConst(i);
175 
176  double sum = 0.0; // do not use Scalar here (Stokhos)
177  for (size_t colID = 0; colID < static_cast<size_t>(aggregate.length); colID++)
178  sum += fineCoordsRandomView(aggregate(colID),j);
179 
180  coarseCoordsView(i,j) = sum / aggregate.length;
181  });
182  }
183  }
184 
185  Set<RCP<doubleMultiVector> >(coarseLevel, "Coordinates", coarseCoords);
186 
187  const ParameterList& pL = GetParameterList();
188  int writeStart = pL.get<int>("write start"), writeEnd = pL.get<int>("write end");
189  if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd) {
190  std::string fileName = "coordinates_before_rebalance_level_" + toString(fineLevel.GetLevelID()) + ".m";
191  Xpetra::IO<double,LO,GO,NO>::Write(fileName, *fineCoords);
192  }
193  if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd) {
194  std::string fileName = "coordinates_before_rebalance_level_" + toString(coarseLevel.GetLevelID()) + ".m";
195  Xpetra::IO<double,LO,GO,NO>::Write(fileName,*coarseCoords);
196  }
197  }
198 
199 } // namespace MueLu
200 
201 #endif // MUELU_COORDINATESTRANSFER_FACTORY_KOKKOS_DEF_HPP
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
One-liner description of what is happening.
Namespace for MueLu classes and methods.