MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_IndexManager_kokkos_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_INDEXMANAGER_KOKKOS_DECL_HPP
11#define MUELU_INDEXMANAGER_KOKKOS_DECL_HPP
12
13#include "MueLu_ConfigDefs.hpp"
14#include "MueLu_Types.hpp"
15
16#include <Tpetra_KokkosCompat_ClassicNodeAPI_Wrapper.hpp>
17
18#include "Teuchos_OrdinalTraits.hpp"
19
20#include "MueLu_BaseClass.hpp"
22
23/*****************************************************************************
24
25****************************************************************************/
26
27namespace MueLu {
28
40template <class LocalOrdinal, class GlobalOrdinal, class Node>
42#undef MUELU_INDEXMANAGER_KOKKOS_SHORT
44
45 public:
46 using execution_space = typename Node::execution_space;
47 using memory_space = typename Node::memory_space;
48 using device_type = Kokkos::Device<execution_space, memory_space>;
49 using intTupleView = typename Kokkos::View<int[3], device_type>;
50 using LOTupleView = typename Kokkos::View<LO[3], device_type>;
51
52 private:
53 const int meshLayout = UNCOUPLED;
54 int myRank = -1;
59
63
67
68 public:
71
73 IndexManager_kokkos(const int NumDimensions,
74 const int interpolationOrder,
75 const int MyRank,
76 const ArrayView<const LO> LFineNodesPerDir,
77 const ArrayView<const int> CoarseRate);
78
80
82 void setupIM(const int NumDimensions,
83 const int interpolationOrder,
84 const ArrayView<const int> coarseRate,
85 const ArrayView<const LO> LFineNodesPerDir);
86
90
91 int getNumDimensions() const { return numDimensions; }
92
94
95 LO getNumLocalFineNodes() const { return lNumFineNodes; }
96
97 LO getNumCoarseNodes() const { return numCoarseNodes; }
98
99 KOKKOS_INLINE_FUNCTION
101
102 KOKKOS_INLINE_FUNCTION
104
105 KOKKOS_INLINE_FUNCTION
107
108 KOKKOS_INLINE_FUNCTION
110
111 Array<LO> getCoarseNodesPerDirArray() const;
112
113 KOKKOS_INLINE_FUNCTION
114 void getFineLID2FineTuple(const LO myLID, LO (&tuple)[3]) const {
115 LO tmp;
116 tuple[2] = myLID / (lFineNodesPerDir(1) * lFineNodesPerDir(0));
117 tmp = myLID % (lFineNodesPerDir(1) * lFineNodesPerDir(0));
118 tuple[1] = tmp / lFineNodesPerDir(0);
119 tuple[0] = tmp % lFineNodesPerDir(0);
120 } // getFineNodeLocalTuple
121
122 KOKKOS_INLINE_FUNCTION
123 void getFineTuple2FineLID(const LO tuple[3], LO& myLID) const {
124 myLID = tuple[2] * lNumFineNodes10 + tuple[1] * lFineNodesPerDir[0] + tuple[0];
125 } // getFineNodeLID
126
127 KOKKOS_INLINE_FUNCTION
128 void getCoarseLID2CoarseTuple(const LO myLID, LO (&tuple)[3]) const {
129 LO tmp;
130 tuple[2] = myLID / numCoarseNodes10;
131 tmp = myLID % numCoarseNodes10;
132 tuple[1] = tmp / coarseNodesPerDir[0];
133 tuple[0] = tmp % coarseNodesPerDir[0];
134 } // getCoarseNodeLocalTuple
135
136 KOKKOS_INLINE_FUNCTION
137 void getCoarseTuple2CoarseLID(const LO i, const LO j, const LO k, LO& myLID) const {
138 myLID = k * numCoarseNodes10 + j * coarseNodesPerDir[0] + i;
139 } // getCoarseNodeLID
140};
141
142} // namespace MueLu
143
144#define MUELU_INDEXMANAGER_KOKKOS_SHORT
145#endif // MUELU_INDEXMANAGER_KOKKOS_DECL_HPP
Base class for MueLu classes.
Container class for mesh layout and indices calculation.
LO numCoarseNodes10
local number of nodes per 0-1 slice remaining after coarsening.
LOTupleView coarseNodesPerDir
local number of nodes per direction remaing after coarsening.
KOKKOS_INLINE_FUNCTION LOTupleView getCoarseNodesPerDir() const
intTupleView endRate
adapted coarsening rate at the edge of the mesh in each direction.
LO lNumFineNodes10
local number of nodes per 0-1 slice.
typename Kokkos::View< LO[3], device_type > LOTupleView
Kokkos::Device< execution_space, memory_space > device_type
KOKKOS_INLINE_FUNCTION void getCoarseLID2CoarseTuple(const LO myLID, LO(&tuple)[3]) const
KOKKOS_INLINE_FUNCTION void getCoarseTuple2CoarseLID(const LO i, const LO j, const LO k, LO &myLID) const
KOKKOS_INLINE_FUNCTION LOTupleView getLocalFineNodesPerDir() const
typename Node::execution_space execution_space
KOKKOS_INLINE_FUNCTION void getFineLID2FineTuple(const LO myLID, LO(&tuple)[3]) const
LOTupleView lFineNodesPerDir
local number of nodes per direction.
LO numCoarseNodes
local number of nodes remaining after coarsening.
KOKKOS_INLINE_FUNCTION intTupleView getCoarseningRates() const
IndexManager_kokkos()=default
Default constructor, return empty object.
intTupleView coarseRate
coarsening rate in each direction
KOKKOS_INLINE_FUNCTION intTupleView getCoarseningEndRates() const
int numDimensions
Number of spacial dimensions in the problem.
KOKKOS_INLINE_FUNCTION void getFineTuple2FineLID(const LO tuple[3], LO &myLID) 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< int[3], device_type > intTupleView
Namespace for MueLu classes and methods.