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