10 #ifndef MUELU_INDEXMANAGER_DEF_KOKKOS_HPP
11 #define MUELU_INDEXMANAGER_DEF_KOKKOS_HPP
17 #include <Xpetra_MapFactory.hpp>
30 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
33 const int interpolationOrder,
38 , coarseRate(
"coarsening rate")
40 , lFineNodesPerDir(
"lFineNodesPerDir")
41 , coarseNodesPerDir(
"lFineNodesPerDir") {
43 if (
const char* dbg = std::getenv(
"MUELU_INDEXMANAGER_DEBUG")) {
44 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
50 setupIM(NumDimensions, interpolationOrder, CoarseRate, LFineNodesPerDir);
52 *out <<
"Done setting up the IndexManager" << std::endl;
56 *out <<
"Computed Mesh Parameters" << std::endl;
60 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
62 setupIM(
const int NumDimensions,
const int interpolationOrder,
64 numDimensions = NumDimensions;
65 interpolationOrder_ = interpolationOrder;
69 "LFineNodesPerDir has to be of size 3 or of size numDimensions!");
71 typename Kokkos::View<LO[3], device_type>::HostMirror lFineNodesPerDir_h = Kokkos::create_mirror_view(lFineNodesPerDir);
72 Kokkos::deep_copy(lFineNodesPerDir_h, lFineNodesPerDir);
73 typename Kokkos::View<LO[3], device_type>::HostMirror coarseRate_h = Kokkos::create_mirror_view(coarseRate);
74 Kokkos::deep_copy(coarseRate_h, coarseRate);
78 for (
int dim = 0; dim < 3; ++dim) {
79 if (dim < getNumDimensions()) {
80 lFineNodesPerDir_h(dim) = LFineNodesPerDir[dim];
81 if (CoarseRate.
size() == 1) {
82 coarseRate_h(dim) = CoarseRate[0];
83 }
else if (CoarseRate.
size() == getNumDimensions()) {
84 coarseRate_h(dim) = CoarseRate[dim];
87 lFineNodesPerDir_h(dim) = 1;
88 coarseRate_h(dim) = 1;
92 Kokkos::deep_copy(lFineNodesPerDir, lFineNodesPerDir_h);
93 Kokkos::deep_copy(coarseRate, coarseRate_h);
97 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
100 if (
const char* dbg = std::getenv(
"MUELU_INDEXMANAGER_DEBUG")) {
101 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
107 typename Kokkos::View<int[3], device_type>::HostMirror coarseRate_h = Kokkos::create_mirror_view(coarseRate);
108 typename Kokkos::View<int[3], device_type>::HostMirror endRate_h = Kokkos::create_mirror_view(endRate);
110 typename Kokkos::View<LO[3], device_type>::HostMirror lFineNodesPerDir_h = Kokkos::create_mirror_view(lFineNodesPerDir);
111 typename Kokkos::View<LO[3], device_type>::HostMirror coarseNodesPerDir_h = Kokkos::create_mirror_view(coarseNodesPerDir);
112 Kokkos::deep_copy(lFineNodesPerDir_h, lFineNodesPerDir);
113 Kokkos::deep_copy(coarseRate_h, coarseRate);
115 lNumFineNodes10 = lFineNodesPerDir_h(1) * lFineNodesPerDir_h(0);
116 lNumFineNodes = lFineNodesPerDir_h(2) * lNumFineNodes10;
117 for (
int dim = 0; dim < 3; ++dim) {
118 if (dim < numDimensions) {
119 endRate_h(dim) = (lFineNodesPerDir_h(dim) - 1) % coarseRate_h(dim);
120 if (endRate_h(dim) == 0) {
121 endRate_h(dim) = coarseRate_h(dim);
129 *out <<
"lFineNodesPerDir: {" << lFineNodesPerDir_h(0) <<
", " << lFineNodesPerDir_h(1) <<
", "
130 << lFineNodesPerDir_h(2) <<
"}" << std::endl;
131 *out <<
"endRate: {" << endRate_h(0) <<
", " << endRate_h(1) <<
", "
132 << endRate_h(2) <<
"}" << std::endl;
147 for (
int dim = 0; dim < 3; ++dim) {
148 if (dim < numDimensions) {
152 coarseNodesPerDir_h(dim) = (lFineNodesPerDir_h(dim) - endRate_h(dim) - 1) / coarseRate_h(dim) + 2;
156 coarseNodesPerDir_h(dim) = 1;
161 if (lFineNodesPerDir_h(dim) < 1) {
162 coarseNodesPerDir_h(dim) = 0;
167 numCoarseNodes10 = coarseNodesPerDir_h(0) * coarseNodesPerDir_h(1);
168 numCoarseNodes = numCoarseNodes10 * coarseNodesPerDir_h(2);
170 *out <<
"coarseNodesPerDir: {" << coarseNodesPerDir_h(0) <<
", "
171 << coarseNodesPerDir_h(1) <<
", " << coarseNodesPerDir_h(2) <<
"}" << std::endl;
172 *out <<
"numCoarseNodes=" << numCoarseNodes << std::endl;
175 Kokkos::deep_copy(coarseRate, coarseRate_h);
176 Kokkos::deep_copy(endRate, endRate_h);
177 Kokkos::deep_copy(lFineNodesPerDir, lFineNodesPerDir_h);
178 Kokkos::deep_copy(coarseNodesPerDir, coarseNodesPerDir_h);
181 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
184 typename LOTupleView::HostMirror coarseNodesPerDir_h = Kokkos::create_mirror_view(coarseNodesPerDir);
185 Kokkos::deep_copy(coarseNodesPerDir_h, coarseNodesPerDir);
188 for (
int dim = 0; dim < 3; ++dim) {
189 coarseNodesPerDirArray[dim] = coarseNodesPerDir_h(dim);
192 return coarseNodesPerDirArray;
197 #define MUELU_INDEXMANAGER_KOKKOS_SHORT
198 #endif // MUELU_INDEXMANAGER_DEF_KOKKOS_HPP
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.
basic_FancyOStream & setShowProcRank(const bool showProcRank)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void computeMeshParameters()
Array< LO > getCoarseNodesPerDirArray() const
basic_FancyOStream & setShowAllFrontMatter(const bool showAllFrontMatter)
Exception throws to report errors in the internal logical of the program.
IndexManager_kokkos()=default
Default constructor, return empty object.