46 #ifndef MUELU_INDEXMANAGER_KOKKOS_DECL_HPP
47 #define MUELU_INDEXMANAGER_KOKKOS_DECL_HPP
52 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
53 #include <KokkosCompat_ClassicNodeAPI_Wrapper.hpp>
54 #include <Kokkos_View.hpp>
58 #include <Xpetra_Map_fwd.hpp>
59 #include <Xpetra_Vector_fwd.hpp>
60 #include <Xpetra_VectorFactory_fwd.hpp>
82 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
83 class IndexManager_kokkos :
public BaseClass {
84 #undef MUELU_INDEXMANAGER_KOKKOS_SHORT
88 typedef typename Node::execution_space execution_space;
89 typedef typename execution_space::memory_space memory_space;
98 int interpolationOrder_;
99 intTupleView coarseRate;
100 intTupleView endRate;
104 LOTupleView lFineNodesPerDir;
108 LOTupleView coarseNodesPerDir;
113 IndexManager_kokkos() =
default;
116 IndexManager_kokkos(
const int NumDimensions,
117 const int interpolationOrder,
119 const ArrayView<const LO> LFineNodesPerDir,
120 const ArrayView<const int> CoarseRate);
122 virtual ~IndexManager_kokkos() {}
125 void setupIM(
const int NumDimensions,
126 const int interpolationOrder,
127 const ArrayView<const int> coarseRate,
128 const ArrayView<const LO> LFineNodesPerDir);
132 void computeMeshParameters();
134 int getNumDimensions()
const {
return numDimensions;}
136 int getInterpolationOrder()
const {
return interpolationOrder_;}
138 LO getNumLocalFineNodes()
const {
return lNumFineNodes;}
140 LO getNumCoarseNodes()
const {
return numCoarseNodes;}
142 KOKKOS_INLINE_FUNCTION
143 intTupleView getCoarseningRates()
const {
return coarseRate;}
145 KOKKOS_INLINE_FUNCTION
146 intTupleView getCoarseningEndRates()
const {
return endRate;}
148 KOKKOS_INLINE_FUNCTION
149 LOTupleView getLocalFineNodesPerDir()
const {
return lFineNodesPerDir;}
151 KOKKOS_INLINE_FUNCTION
152 LOTupleView getCoarseNodesPerDir()
const {
return coarseNodesPerDir;}
154 Array<LO> getCoarseNodesPerDirArray()
const;
156 KOKKOS_INLINE_FUNCTION
157 void getFineLID2FineTuple(
const LO myLID, LO (&tuple)[3])
const {
159 tuple[2] = myLID / (lFineNodesPerDir(1)*lFineNodesPerDir(0));
160 tmp = myLID % (lFineNodesPerDir(1)*lFineNodesPerDir(0));
161 tuple[1] = tmp / lFineNodesPerDir(0);
162 tuple[0] = tmp % lFineNodesPerDir(0);
165 KOKKOS_INLINE_FUNCTION
166 void getFineTuple2FineLID(
const LO tuple[3], LO& myLID)
const {
167 myLID = tuple[2]*lNumFineNodes10 + tuple[1]*lFineNodesPerDir[0] + tuple[0];
170 KOKKOS_INLINE_FUNCTION
171 void getCoarseLID2CoarseTuple(
const LO myLID, LO (&tuple)[3])
const {
173 tuple[2] = myLID / numCoarseNodes10;
174 tmp = myLID % numCoarseNodes10;
175 tuple[1] = tmp / coarseNodesPerDir[0];
176 tuple[0] = tmp % coarseNodesPerDir[0];
179 KOKKOS_INLINE_FUNCTION
180 void getCoarseTuple2CoarseLID(
const LO i,
const LO j,
const LO k, LO& myLID)
const {
181 myLID = k*numCoarseNodes10 + j*coarseNodesPerDir[0] + i;
188 #define MUELU_INDEXMANAGER_KOKKOS_SHORT
189 #endif // HAVE_MUELU_KOKKOS_REFACTOR
190 #endif // MUELU_INDEXMANAGER_KOKKOS_DECL_HPP