MueLu Version of the Day
Loading...
Searching...
No Matches
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
61namespace 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
86 return validParamList;
87 }
88
89 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
91 static bool isAvailableCoords = false;
92
93 const ParameterList& pL = GetParameterList();
94 if(pL.get<bool>("structured aggregation") == true) {
95 if(pL.get<bool>("aggregation coupled") == true) {
96 Input(fineLevel, "gCoarseNodesPerDim");
97 }
98 Input(fineLevel, "lCoarseNodesPerDim");
99 Input(fineLevel, "numDimensions");
100 } else if(pL.get<bool>("Geometric") == true) {
101 Input(coarseLevel, "coarseCoordinates");
102 Input(coarseLevel, "gCoarseNodesPerDim");
103 Input(coarseLevel, "lCoarseNodesPerDim");
104 } else if(pL.get<bool>("hybrid aggregation") == true) {
105 Input(fineLevel, "aggregationRegionTypeCoarse");
106 Input(fineLevel, "lCoarseNodesPerDim");
107 Input(fineLevel, "numDimensions");
108 if(pL.get<bool>("interface aggregation") == true) {
109 Input(fineLevel, "coarseInterfacesDimensions");
110 Input(fineLevel, "nodeOnCoarseInterface");
111 }
112 } else {
113 if (coarseLevel.GetRequestMode() == Level::REQUEST)
114 isAvailableCoords = coarseLevel.IsAvailable("Coordinates", this);
115
116 if (isAvailableCoords == false) {
117 Input(fineLevel, "Coordinates");
118 Input(fineLevel, "Aggregates");
119 Input(fineLevel, "CoarseMap");
120 }
121 }
122 }
123
124 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
126 FactoryMonitor m(*this, "Build", coarseLevel);
127
128 using xdMV = Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO>;
129
130 GetOStream(Runtime0) << "Transferring coordinates" << std::endl;
131
132 int numDimensions;
133 RCP<xdMV> coarseCoords;
134 RCP<xdMV> fineCoords;
135 Array<GO> gCoarseNodesPerDir;
136 Array<LO> lCoarseNodesPerDir;
137
138 const ParameterList& pL = GetParameterList();
139
140 if(pL.get<bool>("hybrid aggregation") == true) {
141 std::string regionType = Get<std::string>(fineLevel,"aggregationRegionTypeCoarse");
142 numDimensions = Get<int>(fineLevel, "numDimensions");
143 lCoarseNodesPerDir = Get<Array<LO> >(fineLevel, "lCoarseNodesPerDim");
144 Set<std::string>(coarseLevel, "aggregationRegionType", regionType);
145 Set<int> (coarseLevel, "numDimensions", numDimensions);
146 Set<Array<LO> > (coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
147
148 if((pL.get<bool>("interface aggregation") == true) && (regionType == "uncoupled")) {
149 Array<LO> coarseInterfacesDimensions = Get<Array<LO> >(fineLevel, "coarseInterfacesDimensions");
150 Array<LO> nodeOnCoarseInterface = Get<Array<LO> >(fineLevel, "nodeOnCoarseInterface");
151 Set<Array<LO> >(coarseLevel, "interfacesDimensions", coarseInterfacesDimensions);
152 Set<Array<LO> >(coarseLevel, "nodeOnInterface", nodeOnCoarseInterface);
153 }
154
155 } else if(pL.get<bool>("structured aggregation") == true) {
156 if(pL.get<bool>("aggregation coupled") == true) {
157 gCoarseNodesPerDir = Get<Array<GO> >(fineLevel, "gCoarseNodesPerDim");
158 Set<Array<GO> >(coarseLevel, "gNodesPerDim", gCoarseNodesPerDir);
159 }
160 lCoarseNodesPerDir = Get<Array<LO> >(fineLevel, "lCoarseNodesPerDim");
161 Set<Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
162 numDimensions = Get<int>(fineLevel, "numDimensions");
163 Set<int>(coarseLevel, "numDimensions", numDimensions);
164
165 } else if(pL.get<bool>("Geometric") == true) {
166 coarseCoords = Get<RCP<xdMV> >(coarseLevel, "coarseCoordinates");
167 gCoarseNodesPerDir = Get<Array<GO> >(coarseLevel, "gCoarseNodesPerDim");
168 lCoarseNodesPerDir = Get<Array<LO> >(coarseLevel, "lCoarseNodesPerDim");
169 Set<Array<GO> >(coarseLevel, "gNodesPerDim", gCoarseNodesPerDir);
170 Set<Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
171
172 Set<RCP<xdMV> >(coarseLevel, "Coordinates", coarseCoords);
173
174 } else {
175 if (coarseLevel.IsAvailable("Coordinates", this)) {
176 GetOStream(Runtime0) << "Reusing coordinates" << std::endl;
177 return;
178 }
179
180 fineCoords = Get< RCP<xdMV> >(fineLevel, "Coordinates");
181 RCP<const Map> coarseMap = Get< RCP<const Map> > (fineLevel, "CoarseMap");
182
183 // coarseMap is being used to set up the domain map of tentative P, and therefore, the row map of Ac
184 // Therefore, if we amalgamate coarseMap, logical nodes in the coordinates vector would correspond to
185 // logical blocks in the matrix
186
187 ArrayView<const GO> elementAList = coarseMap->getLocalElementList();
188
189 LO blkSize = 1;
190 if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null)
191 blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap)->getFixedBlockSize();
192
193 GO indexBase = coarseMap->getIndexBase();
194 size_t numElements = elementAList.size() / blkSize;
195 Array<GO> elementList(numElements);
196
197 // Amalgamate the map
198 for (LO i = 0; i < Teuchos::as<LO>(numElements); i++)
199 elementList[i] = (elementAList[i*blkSize]-indexBase)/blkSize + indexBase;
200
201 RCP<const Map> uniqueMap = fineCoords->getMap();
202 RCP<const Map> coarseCoordMap = MapFactory ::Build(coarseMap->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), elementList, indexBase, coarseMap->getComm());
203 coarseCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO>::Build(coarseCoordMap, fineCoords->getNumVectors());
204
205
206 RCP<Aggregates> aggregates;
207 bool aggregatesCrossProcessors;
208 aggregates = Get<RCP<Aggregates> >(fineLevel, "Aggregates");
209 aggregatesCrossProcessors = aggregates->AggregatesCrossProcessors();
210
211 // Create overlapped fine coordinates to reduce global communication
212 RCP<xdMV> ghostedCoords = fineCoords;
213 if (aggregatesCrossProcessors) {
214 RCP<const Map> nonUniqueMap = aggregates->GetMap();
215 RCP<const Import> importer = ImportFactory::Build(uniqueMap, nonUniqueMap);
216
217 ghostedCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO>::Build(nonUniqueMap, fineCoords->getNumVectors());
218 ghostedCoords->doImport(*fineCoords, *importer, Xpetra::INSERT);
219 }
220
221 // The good news is that this graph has already been constructed for the
222 // TentativePFactory and was cached in Aggregates. So this is a no-op.
223 auto aggGraph = aggregates->GetGraph();
224 auto numAggs = aggGraph.numRows();
225
226 auto fineCoordsView = ghostedCoords->getDeviceLocalView(Xpetra::Access::ReadOnly);
227 auto coarseCoordsView = coarseCoords->getDeviceLocalView(Xpetra::Access::OverwriteAll);
228
229 // Fill in coarse coordinates
230 {
231 SubFactoryMonitor m2(*this, "AverageCoords", coarseLevel);
232
233 const auto dim = ghostedCoords->getNumVectors();
234
235 typename AppendTrait<decltype(fineCoordsView), Kokkos::RandomAccess>::type fineCoordsRandomView = fineCoordsView;
236 for (size_t j = 0; j < dim; j++) {
237 Kokkos::parallel_for("MueLu:CoordinatesTransferF:Build:coord", Kokkos::RangePolicy<local_ordinal_type, execution_space>(0, numAggs),
238 KOKKOS_LAMBDA(const LO i) {
239 // A row in this graph represents all node ids in the aggregate
240 // Therefore, averaging is very easy
241
242 auto aggregate = aggGraph.rowConst(i);
243
244 typename Teuchos::ScalarTraits<Scalar>::magnitudeType sum = 0.0; // do not use Scalar here (Stokhos)
245 for (size_t colID = 0; colID < static_cast<size_t>(aggregate.length); colID++)
246 sum += fineCoordsRandomView(aggregate(colID),j);
247
248 coarseCoordsView(i,j) = sum / aggregate.length;
249 });
250 }
251 }
252
253 Set<RCP<xdMV> >(coarseLevel, "Coordinates", coarseCoords);
254
255 }
256
257 int writeStart = pL.get<int>("write start"), writeEnd = pL.get<int>("write end");
258 if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd) {
259 std::ostringstream buf;
260 buf << fineLevel.GetLevelID();
261 std::string fileName = "coordinates_before_rebalance_level_" + buf.str() + ".m";
262 Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO>::Write(fileName,*fineCoords);
263 }
264 if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd) {
265 std::ostringstream buf;
266 buf << coarseLevel.GetLevelID();
267 std::string fileName = "coordinates_before_rebalance_level_" + buf.str() + ".m";
268 Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO>::Write(fileName,*coarseCoords);
269 }
270 }
271
272} // namespace MueLu
273
274#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.
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.