MueLu  Version of the Day
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_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< 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 Node>
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 Node>
93  FactoryMonitor m(*this, "Build", coarseLevel);
94 
95  typedef Xpetra::MultiVector<double,LO,GO,NO> xdMV;
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  RCP<Aggregates> aggregates = Get< RCP<Aggregates> > (fineLevel, "Aggregates");
105  RCP<xdMV> fineCoords = Get< RCP<xdMV> >(fineLevel, "Coordinates");
106  RCP<const Map> 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  LO blkSize = 1;
114  if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null)
115  blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap)->getFixedBlockSize();
116 
117  GO indexBase = coarseMap->getIndexBase();
118  size_t numElements = elementAList.size() / blkSize;
119  Array<GO> elementList(numElements);
120 
121  // Amalgamate the map
122  for (LO i = 0; i < Teuchos::as<LO>(numElements); i++)
123  elementList[i] = (elementAList[i*blkSize]-indexBase)/blkSize + indexBase;
124 
125  RCP<const Map> uniqueMap = fineCoords->getMap();
126  RCP<const Map> coarseCoordMap = MapFactory ::Build(coarseMap->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), elementList, indexBase, coarseMap->getComm());
127  RCP<xdMV> coarseCoords = Xpetra::MultiVectorFactory<double,LO,GO,NO>::Build(coarseCoordMap, fineCoords->getNumVectors());
128 
129  // Create overlapped fine coordinates to reduce global communication
130  RCP<xdMV> ghostedCoords = fineCoords;
131  if (aggregates->AggregatesCrossProcessors()) {
132  RCP<const Map> nonUniqueMap = aggregates->GetMap();
133  RCP<const Import> importer = ImportFactory::Build(uniqueMap, nonUniqueMap);
134 
135  ghostedCoords = Xpetra::MultiVectorFactory<double,LO,GO,NO>::Build(nonUniqueMap, fineCoords->getNumVectors());
136  ghostedCoords->doImport(*fineCoords, *importer, Xpetra::INSERT);
137  }
138 
139  // Get some info about aggregates
140  int myPID = uniqueMap->getComm()->getRank();
141  LO numAggs = aggregates->GetNumAggregates();
142  ArrayRCP<LO> aggSizes = aggregates->ComputeAggregateSizes(true,true);
143  const ArrayRCP<const LO> vertex2AggID = aggregates->GetVertex2AggId()->getData(0);
144  const ArrayRCP<const LO> procWinner = aggregates->GetProcWinner()->getData(0);
145 
146  // Fill in coarse coordinates
147  for (size_t j = 0; j < fineCoords->getNumVectors(); j++) {
148  ArrayRCP<const double> fineCoordsData = ghostedCoords->getData(j);
149  ArrayRCP<double> coarseCoordsData = coarseCoords->getDataNonConst(j);
150 
151  for (LO lnode = 0; lnode < vertex2AggID.size(); lnode++)
152  if (procWinner[lnode] == myPID)
153  coarseCoordsData[vertex2AggID[lnode]] += fineCoordsData[lnode];
154 
155  for (LO agg = 0; agg < numAggs; agg++)
156  coarseCoordsData[agg] /= aggSizes[agg];
157  }
158 
159  Set<RCP<xdMV> >(coarseLevel, "Coordinates", coarseCoords);
160 
161  const ParameterList& pL = GetParameterList();
162  int writeStart = pL.get<int>("write start"), writeEnd = pL.get<int>("write end");
163  if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd) {
164  std::ostringstream buf;
165  buf << fineLevel.GetLevelID();
166  std::string fileName = "coordinates_before_rebalance_level_" + buf.str() + ".m";
167  Xpetra::IO<double,LO,GO,NO>::Write(fileName,*fineCoords);
168  }
169  if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd) {
170  std::ostringstream buf;
171  buf << coarseLevel.GetLevelID();
172  std::string fileName = "coordinates_before_rebalance_level_" + buf.str() + ".m";
173  Xpetra::IO<double,LO,GO,NO>::Write(fileName,*coarseCoords);
174  }
175  }
176 
177 } // namespace MueLu
178 
179 #endif // MUELU_COORDINATESTRANSFER_FACTORY_DEF_HPP
size_type size() const
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
RequestMode GetRequestMode() const
Timer to be used in factories. Similar to Monitor but with additional timers.
One-liner description of what is happening.
Namespace for MueLu classes and methods.
void DeclareInput(Level &finelevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const RCP< LOVector > & GetVertex2AggId() const
Returns constant vector that maps local node IDs to local aggregates IDs.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
const RCP< const Map > GetMap() const
returns (overlapping) map of aggregate/node distribution
LO GetNumAggregates() const
returns the number of aggregates of the current processor. Note: could/should be renamed to GetNumLoc...
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need&#39;s value has been saved.
const RCP< LOVector > & GetProcWinner() const
Returns constant vector that maps local node IDs to owning processor IDs.
void AggregatesCrossProcessors(const bool &flag)
Record whether aggregates include DOFs from other processes.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Teuchos::ArrayRCP< LO > ComputeAggregateSizes(bool forceRecompute=true, bool cacheSizes=false) const
Compute sizes of aggregates.