MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_UncoupledIndexManager_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 // Ray Tuminaro (rstumin@sandia.gov)
41 // Luc Berger-Vergoat (lberge@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_UNCOUPLEDINDEXMANAGER_DEF_HPP_
47 #define MUELU_UNCOUPLEDINDEXMANAGER_DEF_HPP_
48 
49 #include <Xpetra_MapFactory.hpp>
52 
53 namespace MueLu {
54 
55 template <class LocalOrdinal, class GlobalOrdinal, class Node>
57  UncoupledIndexManager(const RCP<const Teuchos::Comm<int> > comm, const bool coupled,
58  const int NumDimensions, const int interpolationOrder,
59  const int MyRank, const int NumRanks,
60  const Array<GO> GFineNodesPerDir, const Array<LO> LFineNodesPerDir,
61  const Array<LO> CoarseRate, const bool singleCoarsePoint)
62  : IndexManager(comm, coupled, singleCoarsePoint, NumDimensions, interpolationOrder,
63  Array<GO>(3, -1), LFineNodesPerDir)
64  , myRank(MyRank)
65  , numRanks(NumRanks) {
66  // Load coarse rate, being careful about formating
67  for (int dim = 0; dim < 3; ++dim) {
68  if (dim < this->numDimensions) {
69  if (CoarseRate.size() == 1) {
70  this->coarseRate[dim] = CoarseRate[0];
71  } else if (CoarseRate.size() == this->numDimensions) {
72  this->coarseRate[dim] = CoarseRate[dim];
73  }
74  } else {
75  this->coarseRate[dim] = 1;
76  }
77  }
78 
79  this->computeMeshParameters();
82 } // Constructor
83 
84 template <class LocalOrdinal, class GlobalOrdinal, class Node>
87  GO input[1] = {as<GO>(this->lNumCoarseNodes)}, output[1] = {0};
88  Teuchos::reduceAll(*(this->comm_), Teuchos::REDUCE_SUM, 1, input, output);
89  this->gNumCoarseNodes = output[0];
90 } // computeGlobalCoarseParameters
91 
92 template <class LocalOrdinal, class GlobalOrdinal, class Node>
95  Array<LO>& ghostedNodeCoarseLIDs,
96  Array<int>& ghostedNodeCoarsePIDs,
97  Array<GO>& /* ghostedNodeCoarseGIDs */) const {
98  // First we allocate memory for the outputs
99  ghostedNodeCoarseLIDs.resize(this->getNumLocalGhostedNodes());
100  ghostedNodeCoarsePIDs.resize(this->getNumLocalGhostedNodes());
101  // In the uncoupled case the data required is trivial to provide!
102  for (LO idx = 0; idx < this->getNumLocalGhostedNodes(); ++idx) {
103  ghostedNodeCoarseLIDs[idx] = idx;
104  ghostedNodeCoarsePIDs[idx] = myRank;
105  }
106 } // getGhostedNodesData
107 
108 template <class LocalOrdinal, class GlobalOrdinal, class Node>
110  getCoarseNodesData(const RCP<const Map> fineCoordinatesMap,
111  Array<GO>& coarseNodeCoarseGIDs,
112  Array<GO>& coarseNodeFineGIDs) const {
113  // Allocate sufficient amount of storage in output arrays
114  coarseNodeCoarseGIDs.resize(this->getNumLocalCoarseNodes());
115  coarseNodeFineGIDs.resize(this->getNumLocalCoarseNodes());
116 
117  // Load all the GIDs on the fine mesh
118  ArrayView<const GO> fineNodeGIDs = fineCoordinatesMap->getLocalElementList();
119 
120  // Extract the fine LIDs of the coarse nodes and store the corresponding GIDs
121  LO fineLID;
122  for (LO coarseLID = 0; coarseLID < this->getNumLocalCoarseNodes(); ++coarseLID) {
123  Array<LO> coarseIndices(3), fineIndices(3);
124  this->getCoarseNodeLocalTuple(coarseLID,
125  coarseIndices[0],
126  coarseIndices[1],
127  coarseIndices[2]);
128  for (int dim = 0; dim < 3; ++dim) {
129  if (coarseIndices[dim] == this->lCoarseNodesPerDir[dim] - 1) {
130  if (this->lCoarseNodesPerDir[dim] == 1) {
131  fineIndices[dim] = 0;
132  } else {
133  fineIndices[dim] = this->lFineNodesPerDir[dim] - 1;
134  }
135  } else {
136  fineIndices[dim] = coarseIndices[dim] * this->coarseRate[dim];
137  }
138  }
139 
140  fineLID = fineIndices[2] * this->lNumFineNodes10 + fineIndices[1] * this->lFineNodesPerDir[0] + fineIndices[0];
141  coarseNodeFineGIDs[coarseLID] = fineNodeGIDs[fineLID];
142  }
143 } // getCoarseNodesData
144 
145 template <class LocalOrdinal, class GlobalOrdinal, class Node>
146 std::vector<std::vector<GlobalOrdinal> > UncoupledIndexManager<LocalOrdinal, GlobalOrdinal, Node>::
148  std::vector<std::vector<GO> > coarseMeshData;
149  return coarseMeshData;
150 }
151 
152 template <class LocalOrdinal, class GlobalOrdinal, class Node>
154  getFineNodeGlobalTuple(const GO /* myGID */, GO& /* i */, GO& /* j */, GO& /* k */) const {
155 }
156 
157 template <class LocalOrdinal, class GlobalOrdinal, class Node>
159  getFineNodeLocalTuple(const LO myLID, LO& i, LO& j, LO& k) const {
160  LO tmp;
161  k = myLID / this->lNumFineNodes10;
162  tmp = myLID % this->lNumFineNodes10;
163  j = tmp / this->lFineNodesPerDir[0];
164  i = tmp % this->lFineNodesPerDir[0];
165 } // getFineNodeLocalTuple
166 
167 template <class LocalOrdinal, class GlobalOrdinal, class Node>
169  getFineNodeGhostedTuple(const LO myLID, LO& i, LO& j, LO& k) const {
170  LO tmp;
171  k = myLID / this->lNumFineNodes10;
172  tmp = myLID % this->lNumFineNodes10;
173  j = tmp / this->lFineNodesPerDir[0];
174  i = tmp % this->lFineNodesPerDir[0];
175 
176  k += this->offsets[2];
177  j += this->offsets[1];
178  i += this->offsets[0];
179 } // getFineNodeGhostedTuple
180 
181 template <class LocalOrdinal, class GlobalOrdinal, class Node>
183  getFineNodeGID(const GO /* i */, const GO /* j */, const GO /* k */, GO& /* myGID */) const {
184 }
185 
186 template <class LocalOrdinal, class GlobalOrdinal, class Node>
188  getFineNodeLID(const LO /* i */, const LO /* j */, const LO /* k */, LO& /* myLID */) const {
189 }
190 
191 template <class LocalOrdinal, class GlobalOrdinal, class Node>
193  getCoarseNodeGlobalTuple(const GO /* myGID */, GO& /* i */, GO& /* j */, GO& /* k */) const {
194 }
195 
196 template <class LocalOrdinal, class GlobalOrdinal, class Node>
198  getCoarseNodeLocalTuple(const LO myLID, LO& i, LO& j, LO& k) const {
199  LO tmp;
200  k = myLID / this->lNumCoarseNodes10;
201  tmp = myLID % this->lNumCoarseNodes10;
202  j = tmp / this->lCoarseNodesPerDir[0];
203  i = tmp % this->lCoarseNodesPerDir[0];
204 } // getCoarseNodeLocalTuple
205 
206 template <class LocalOrdinal, class GlobalOrdinal, class Node>
208  getCoarseNodeGID(const GO /* i */, const GO /* j */, const GO /* k */, GO& /* myGID */) const {
209 }
210 
211 template <class LocalOrdinal, class GlobalOrdinal, class Node>
213  getCoarseNodeLID(const LO /* i */, const LO /* j */, const LO /* k */, LO& /* myLID */) const {
214 }
215 
216 template <class LocalOrdinal, class GlobalOrdinal, class Node>
218  getCoarseNodeGhostedLID(const LO i, const LO j, const LO k, LO& myLID) const {
219  myLID = k * this->numGhostedNodes10 + j * this->ghostedNodesPerDir[0] + i;
220 } // getCoarseNodeGhostedLID
221 
222 template <class LocalOrdinal, class GlobalOrdinal, class Node>
224  getCoarseNodeFineLID(const LO /* i */, const LO /* j */, const LO /* k */, LO& /* myLID */) const {
225 }
226 
227 template <class LocalOrdinal, class GlobalOrdinal, class Node>
229  getGhostedNodeFineLID(const LO /* i */, const LO /* j */, const LO /* k */, LO& /* myLID */) const {
230 }
231 
232 template <class LocalOrdinal, class GlobalOrdinal, class Node>
234  getGhostedNodeCoarseLID(const LO /* i */, const LO /* j */, const LO /* k */, LO& /* myLID */) const {
235 }
236 
237 } // namespace MueLu
238 
239 #endif /* MUELU_UNCOUPLEDINDEXMANAGER_DEF_HPP_ */
void getGhostedNodeCoarseLID(const LO i, const LO j, const LO k, LO &myLID) const
void getFineNodeLID(const LO i, const LO j, const LO k, LO &myLID) const
GlobalOrdinal GO
LocalOrdinal LO
void getCoarseNodeGID(const GO i, const GO j, const GO k, GO &myGID) const
Array< int > coarseRate
coarsening rate in each direction
void getFineNodeLocalTuple(const LO myLID, LO &i, LO &j, LO &k) const
void getFineNodeGhostedTuple(const LO myLID, LO &i, LO &j, LO &k) const
const int numDimensions
Number of spacial dimensions in the problem.
void resize(size_type new_size, const value_type &x=value_type())
void getFineNodeGID(const GO i, const GO j, const GO k, GO &myGID) const
void getCoarseNodeLocalTuple(const LO myLID, LO &i, LO &j, LO &k) const
void getCoarseNodeFineLID(const LO i, const LO j, const LO k, LO &myLID) const
void getCoarseNodeGhostedLID(const LO i, const LO j, const LO k, LO &myLID) const
size_type size() const
GO gNumCoarseNodes10
global number of nodes per 0-1 slice remaining after coarsening.
void getCoarseNodeGlobalTuple(const GO myGID, GO &i, GO &j, GO &k) const
void getFineNodeGlobalTuple(const GO myGID, GO &i, GO &j, GO &k) const
std::vector< std::vector< GO > > getCoarseMeshData() const
void getGhostedNodeFineLID(const LO i, const LO j, const LO k, LO &myLID) const
GO gNumCoarseNodes
global number of nodes remaining after coarsening.
void getCoarseNodesData(const RCP< const Map > fineCoordinatesMap, Array< GO > &coarseNodeCoarseGIDs, Array< GO > &coarseNodeFineGIDs) const
Container class for mesh layout and indices calculation.
void getGhostedNodesData(const RCP< const Map > fineMap, Array< LO > &ghostedNodeCoarseLIDs, Array< int > &ghostedNodeCoarsePIDs, Array< GO > &ghostedNodeCoarseGIDs) const
void getCoarseNodeLID(const LO i, const LO j, const LO k, LO &myLID) const