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"
51 #include "Xpetra_MapFactory.hpp"
52 #include "Xpetra_IO.hpp"
53 
54 #include "MueLu_CoarseMapFactory.hpp"
55 #include "MueLu_Aggregates.hpp"
57 //#include "MueLu_Utilities.hpp"
58 
59 #include "MueLu_Level.hpp"
60 #include "MueLu_Monitor.hpp"
61 
62 namespace MueLu {
63 
64  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
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<bool> ("structured aggregation", false, "Flag specifying that the geometric data is transferred for StructuredAggregationFactory");
72  validParamList->set<bool> ("aggregation coupled", false, "Flag specifying if the aggregation algorithm was used in coupled mode.");
73  validParamList->set<bool> ("Geometric", false, "Flag specifying that the coordinates are transferred for GeneralGeometricPFactory");
74  validParamList->set<RCP<const FactoryBase> >("coarseCoordinates", Teuchos::null, "Factory for coarse coordinates generation");
75  validParamList->set<RCP<const FactoryBase> >("gCoarseNodesPerDim", Teuchos::null, "Factory providing the global number of nodes per spatial dimensions of the mesh");
76  validParamList->set<RCP<const FactoryBase> >("lCoarseNodesPerDim", Teuchos::null, "Factory providing the local number of nodes per spatial dimensions of the mesh");
77  validParamList->set<RCP<const FactoryBase> >("numDimensions" , Teuchos::null, "Factory providing the number of spatial dimensions of the mesh");
78  validParamList->set<int> ("write start", -1, "first level at which coordinates should be written to file");
79  validParamList->set<int> ("write end", -1, "last level at which coordinates should be written to file");
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> ("hybrid aggregation", false, "Flag specifying that hybrid aggregation data is transfered for HybridAggregationFactory");
82 
83 
84  return validParamList;
85  }
86 
87  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
89  static bool isAvailableCoords = false;
90 
91  const ParameterList& pL = GetParameterList();
92  if(pL.get<bool>("structured aggregation") == true) {
93  if(pL.get<bool>("aggregation coupled") == true) {
94  Input(fineLevel, "gCoarseNodesPerDim");
95  }
96  Input(fineLevel, "lCoarseNodesPerDim");
97  Input(fineLevel, "numDimensions");
98  } else if(pL.get<bool>("Geometric") == true) {
99  Input(coarseLevel, "coarseCoordinates");
100  Input(coarseLevel, "gCoarseNodesPerDim");
101  Input(coarseLevel, "lCoarseNodesPerDim");
102  } else {
103  if (coarseLevel.GetRequestMode() == Level::REQUEST)
104  isAvailableCoords = coarseLevel.IsAvailable("Coordinates", this);
105 
106  if (isAvailableCoords == false) {
107  Input(fineLevel, "Coordinates");
108  Input(fineLevel, "Aggregates");
109  Input(fineLevel, "CoarseMap");
110  }
111  }
112  if(pL.get<bool>("hybrid aggregation") == true) {
113  Input(fineLevel,"aggregationRegionTypeCoarse");
114  Input(fineLevel, "lCoarseNodesPerDim");
115  }
116  }
117 
118  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
120  FactoryMonitor m(*this, "Build", coarseLevel);
121 
123 
124  GetOStream(Runtime0) << "Transferring coordinates" << std::endl;
125 
126  int numDimensions;
127  RCP<xdMV> coarseCoords;
128  RCP<xdMV> fineCoords;
129  Array<GO> gCoarseNodesPerDir;
130  Array<LO> lCoarseNodesPerDir;
131 
132  const ParameterList& pL = GetParameterList();
133 
134  if(pL.get<bool>("hybrid aggregation") == true) {
135  std::string regionType = Get<std::string>(fineLevel,"aggregationRegionTypeCoarse");
136  lCoarseNodesPerDir = Get<Array<LO> >(fineLevel, "lCoarseNodesPerDim");
137  Set<std::string>(coarseLevel, "aggregationRegionType", regionType);
138  Set< Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
139  }
140 
141  if(pL.get<bool>("structured aggregation") == true) {
142  if(pL.get<bool>("aggregation coupled") == true) {
143  gCoarseNodesPerDir = Get<Array<GO> >(fineLevel, "gCoarseNodesPerDim");
144  Set<Array<GO> >(coarseLevel, "gNodesPerDim", gCoarseNodesPerDir);
145  }
146  lCoarseNodesPerDir = Get<Array<LO> >(fineLevel, "lCoarseNodesPerDim");
147  Set<Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
148  numDimensions = Get<int>(fineLevel, "numDimensions");
149  Set<int>(coarseLevel, "numDimensions", numDimensions);
150  } else if(pL.get<bool>("Geometric") == true) {
151  coarseCoords = Get<RCP<xdMV> >(coarseLevel, "coarseCoordinates");
152  gCoarseNodesPerDir = Get<Array<GO> >(coarseLevel, "gCoarseNodesPerDim");
153  lCoarseNodesPerDir = Get<Array<LO> >(coarseLevel, "lCoarseNodesPerDim");
154  Set<Array<GO> >(coarseLevel, "gNodesPerDim", gCoarseNodesPerDir);
155  Set<Array<LO> >(coarseLevel, "lNodesPerDim", lCoarseNodesPerDir);
156 
157  Set<RCP<xdMV> >(coarseLevel, "Coordinates", coarseCoords);
158 
159  } else {
160  if (coarseLevel.IsAvailable("Coordinates", this)) {
161  GetOStream(Runtime0) << "Reusing coordinates" << std::endl;
162  return;
163  }
164 
165  RCP<Aggregates> aggregates = Get< RCP<Aggregates> > (fineLevel, "Aggregates");
166  fineCoords = Get< RCP<xdMV> >(fineLevel, "Coordinates");
167  RCP<const Map> coarseMap = Get< RCP<const Map> > (fineLevel, "CoarseMap");
168 
169  // coarseMap is being used to set up the domain map of tentative P, and therefore, the row map of Ac
170  // Therefore, if we amalgamate coarseMap, logical nodes in the coordinates vector would correspond to
171  // logical blocks in the matrix
172 
173  ArrayView<const GO> elementAList = coarseMap->getNodeElementList();
174 
175  LO blkSize = 1;
176  if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null)
177  blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap)->getFixedBlockSize();
178 
179  GO indexBase = coarseMap->getIndexBase();
180  size_t numElements = elementAList.size() / blkSize;
181  Array<GO> elementList(numElements);
182 
183  // Amalgamate the map
184  for (LO i = 0; i < Teuchos::as<LO>(numElements); i++)
185  elementList[i] = (elementAList[i*blkSize]-indexBase)/blkSize + indexBase;
186 
187  RCP<const Map> uniqueMap = fineCoords->getMap();
188  RCP<const Map> coarseCoordMap = MapFactory ::Build(coarseMap->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), elementList, indexBase, coarseMap->getComm());
189  coarseCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LO,GO,NO>::Build(coarseCoordMap, fineCoords->getNumVectors());
190 
191  // Create overlapped fine coordinates to reduce global communication
192  RCP<xdMV> ghostedCoords = fineCoords;
193  if (aggregates->AggregatesCrossProcessors()) {
194  RCP<const Map> nonUniqueMap = aggregates->GetMap();
195  RCP<const Import> importer = ImportFactory::Build(uniqueMap, nonUniqueMap);
196 
197  ghostedCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LO,GO,NO>::Build(nonUniqueMap, fineCoords->getNumVectors());
198  ghostedCoords->doImport(*fineCoords, *importer, Xpetra::INSERT);
199  }
200 
201  // Get some info about aggregates
202  int myPID = uniqueMap->getComm()->getRank();
203  LO numAggs = aggregates->GetNumAggregates();
204  ArrayRCP<LO> aggSizes = aggregates->ComputeAggregateSizes();
205  const ArrayRCP<const LO> vertex2AggID = aggregates->GetVertex2AggId()->getData(0);
206  const ArrayRCP<const LO> procWinner = aggregates->GetProcWinner()->getData(0);
207 
208  // Fill in coarse coordinates
209  for (size_t j = 0; j < fineCoords->getNumVectors(); j++) {
210  ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType> fineCoordsData = ghostedCoords->getData(j);
211  ArrayRCP<typename Teuchos::ScalarTraits<Scalar>::coordinateType> coarseCoordsData = coarseCoords->getDataNonConst(j);
212 
213  for (LO lnode = 0; lnode < vertex2AggID.size(); lnode++) {
214  if (procWinner[lnode] == myPID &&
215  lnode < vertex2AggID.size() &&
216  lnode < fineCoordsData.size() && // TAW do not access off-processor coordinates
217  vertex2AggID[lnode] < coarseCoordsData.size() &&
218  Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<Scalar>::coordinateType>::isnaninf(fineCoordsData[lnode]) == false) {
219  coarseCoordsData[vertex2AggID[lnode]] += fineCoordsData[lnode];
220  }
221  }
222  for (LO agg = 0; agg < numAggs; agg++) {
223  coarseCoordsData[agg] /= aggSizes[agg];
224  }
225  }
226 
227  Set<RCP<xdMV> >(coarseLevel, "Coordinates", coarseCoords);
228  } // if pL.get<bool>("Geometric") == true
229 
230  int writeStart = pL.get<int>("write start"), writeEnd = pL.get<int>("write end");
231  if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd) {
232  std::ostringstream buf;
233  buf << fineLevel.GetLevelID();
234  std::string fileName = "coordinates_before_rebalance_level_" + buf.str() + ".m";
235  Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LO,GO,NO>::Write(fileName,*fineCoords);
236  }
237  if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd) {
238  std::ostringstream buf;
239  buf << coarseLevel.GetLevelID();
240  std::string fileName = "coordinates_before_rebalance_level_" + buf.str() + ".m";
241  Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LO,GO,NO>::Write(fileName,*coarseCoords);
242  }
243  }
244 
245 } // namespace MueLu
246 
247 #endif // MUELU_COORDINATESTRANSFER_FACTORY_DEF_HPP
const RCP< LOVector > & GetProcWinner() const
Returns constant vector that maps local node IDs to owning processor IDs.
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.
const RCP< const Map > GetMap() const
returns (overlapping) map of aggregate/node distribution
LocalOrdinal LO
One-liner description of what is happening.
size_type size() const
Teuchos::ArrayRCP< LO > ComputeAggregateSizes(bool forceRecompute=false) const
Compute sizes of aggregates.
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
LO GetNumAggregates() const
returns the number of aggregates of the current processor. Note: could/should be renamed to GetNumLoc...
const RCP< LOMultiVector > & GetVertex2AggId() const
Returns constant vector that maps local node IDs to local aggregates IDs.
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.
void AggregatesCrossProcessors(const bool &flag)
Record whether aggregates include DOFs from other processes.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need&#39;s value has been saved.