MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_IndexManager_kokkos_decl.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-Vergiat (lberge@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_INDEXMANAGER_KOKKOS_DECL_HPP
47 #define MUELU_INDEXMANAGER_KOKKOS_DECL_HPP
48 
49 #include "MueLu_ConfigDefs.hpp"
50 #include "MueLu_Types.hpp"
51 
52 #include <Tpetra_KokkosCompat_ClassicNodeAPI_Wrapper.hpp>
53 
55 
56 #include "MueLu_BaseClass.hpp"
58 
59 /*****************************************************************************
60 
61 ****************************************************************************/
62 
63 namespace MueLu {
64 
76 template <class LocalOrdinal, class GlobalOrdinal, class Node>
78 #undef MUELU_INDEXMANAGER_KOKKOS_SHORT
80 
81  public:
82  using execution_space = typename Node::execution_space;
83  using memory_space = typename Node::memory_space;
84  using device_type = Kokkos::Device<execution_space, memory_space>;
85  using intTupleView = typename Kokkos::View<int[3], device_type>;
86  using LOTupleView = typename Kokkos::View<LO[3], device_type>;
87 
88  private:
89  const int meshLayout = UNCOUPLED;
90  int myRank = -1;
95 
99 
103 
104  public:
106  IndexManager_kokkos() = default;
107 
109  IndexManager_kokkos(const int NumDimensions,
110  const int interpolationOrder,
111  const int MyRank,
112  const ArrayView<const LO> LFineNodesPerDir,
113  const ArrayView<const int> CoarseRate);
114 
115  virtual ~IndexManager_kokkos() {}
116 
118  void setupIM(const int NumDimensions,
119  const int interpolationOrder,
121  const ArrayView<const LO> LFineNodesPerDir);
122 
125  void computeMeshParameters();
126 
127  int getNumDimensions() const { return numDimensions; }
128 
130 
132 
133  LO getNumCoarseNodes() const { return numCoarseNodes; }
134 
135  KOKKOS_INLINE_FUNCTION
137 
138  KOKKOS_INLINE_FUNCTION
140 
141  KOKKOS_INLINE_FUNCTION
143 
144  KOKKOS_INLINE_FUNCTION
146 
148 
149  KOKKOS_INLINE_FUNCTION
150  void getFineLID2FineTuple(const LO myLID, LO (&tuple)[3]) const {
151  LO tmp;
152  tuple[2] = myLID / (lFineNodesPerDir(1) * lFineNodesPerDir(0));
153  tmp = myLID % (lFineNodesPerDir(1) * lFineNodesPerDir(0));
154  tuple[1] = tmp / lFineNodesPerDir(0);
155  tuple[0] = tmp % lFineNodesPerDir(0);
156  } // getFineNodeLocalTuple
157 
158  KOKKOS_INLINE_FUNCTION
159  void getFineTuple2FineLID(const LO tuple[3], LO& myLID) const {
160  myLID = tuple[2] * lNumFineNodes10 + tuple[1] * lFineNodesPerDir[0] + tuple[0];
161  } // getFineNodeLID
162 
163  KOKKOS_INLINE_FUNCTION
164  void getCoarseLID2CoarseTuple(const LO myLID, LO (&tuple)[3]) const {
165  LO tmp;
166  tuple[2] = myLID / numCoarseNodes10;
167  tmp = myLID % numCoarseNodes10;
168  tuple[1] = tmp / coarseNodesPerDir[0];
169  tuple[0] = tmp % coarseNodesPerDir[0];
170  } // getCoarseNodeLocalTuple
171 
172  KOKKOS_INLINE_FUNCTION
173  void getCoarseTuple2CoarseLID(const LO i, const LO j, const LO k, LO& myLID) const {
174  myLID = k * numCoarseNodes10 + j * coarseNodesPerDir[0] + i;
175  } // getCoarseNodeLID
176 };
177 
178 } // namespace MueLu
179 
180 #define MUELU_INDEXMANAGER_KOKKOS_SHORT
181 #endif // MUELU_INDEXMANAGER_KOKKOS_DECL_HPP
LO lNumFineNodes10
local number of nodes per 0-1 slice.
KOKKOS_INLINE_FUNCTION void getFineLID2FineTuple(const LO myLID, LO(&tuple)[3]) const
int interpolationOrder_
Interpolation order used by grid transfer operators using these aggregates.
void setupIM(const int NumDimensions, const int interpolationOrder, const ArrayView< const int > coarseRate, const ArrayView< const LO > LFineNodesPerDir)
Common setup pattern used for all the different types of undelying mesh.
typename Kokkos::View< LO[3], device_type > LOTupleView
LOTupleView lFineNodesPerDir
local number of nodes per direction.
typename Node::memory_space memory_space
KOKKOS_INLINE_FUNCTION intTupleView getCoarseningEndRates() const
typename Node::execution_space execution_space
KOKKOS_INLINE_FUNCTION intTupleView getCoarseningRates() const
KOKKOS_INLINE_FUNCTION void getCoarseTuple2CoarseLID(const LO i, const LO j, const LO k, LO &myLID) const
LocalOrdinal LO
KOKKOS_INLINE_FUNCTION void getCoarseLID2CoarseTuple(const LO myLID, LO(&tuple)[3]) const
Kokkos::Device< execution_space, memory_space > device_type
LO numCoarseNodes10
local number of nodes per 0-1 slice remaining after coarsening.
LO numCoarseNodes
local number of nodes remaining after coarsening.
Base class for MueLu classes.
KOKKOS_INLINE_FUNCTION void getFineTuple2FineLID(const LO tuple[3], LO &myLID) const
int numDimensions
Number of spacial dimensions in the problem.
intTupleView coarseRate
coarsening rate in each direction
LOTupleView coarseNodesPerDir
local number of nodes per direction remaing after coarsening.
Container class for mesh layout and indices calculation.
typename Kokkos::View< int[3], device_type > intTupleView
KOKKOS_INLINE_FUNCTION LOTupleView getCoarseNodesPerDir() const
KOKKOS_INLINE_FUNCTION LOTupleView getLocalFineNodesPerDir() const
IndexManager_kokkos()=default
Default constructor, return empty object.
intTupleView endRate
adapted coarsening rate at the edge of the mesh in each direction.