MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_IndexManager_kokkos_def.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_DEF_KOKKOS_HPP
11#define MUELU_INDEXMANAGER_DEF_KOKKOS_HPP
12
13#include <utility>
14
15#include "Teuchos_OrdinalTraits.hpp"
16
17#include <Xpetra_MapFactory.hpp>
18
19#include "MueLu_ConfigDefs.hpp"
20#include "MueLu_Types.hpp"
21#include "MueLu_Exceptions.hpp"
23
24/*****************************************************************************
25
26****************************************************************************/
27
28namespace MueLu {
29
30template <class LocalOrdinal, class GlobalOrdinal, class Node>
32 IndexManager_kokkos(const int NumDimensions,
33 const int interpolationOrder,
34 const int MyRank,
35 const ArrayView<const LO> LFineNodesPerDir,
36 const ArrayView<const int> CoarseRate)
37 : myRank(MyRank)
38 , coarseRate("coarsening rate")
39 , endRate("endRate")
40 , lFineNodesPerDir("lFineNodesPerDir")
41 , coarseNodesPerDir("lFineNodesPerDir") {
42 RCP<Teuchos::FancyOStream> out;
43 if (const char* dbg = std::getenv("MUELU_INDEXMANAGER_DEBUG")) {
44 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
45 out->setShowAllFrontMatter(false).setShowProcRank(true);
46 } else {
47 out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
48 }
49
50 setupIM(NumDimensions, interpolationOrder, CoarseRate, LFineNodesPerDir);
51
52 *out << "Done setting up the IndexManager" << std::endl;
53
55
56 *out << "Computed Mesh Parameters" << std::endl;
57
58} // IndexManager_kokkos Constructor
59
60template <class LocalOrdinal, class GlobalOrdinal, class Node>
62 setupIM(const int NumDimensions, const int interpolationOrder,
63 const ArrayView<const int> CoarseRate, const ArrayView<const LO> LFineNodesPerDir) {
64 numDimensions = NumDimensions;
65 interpolationOrder_ = interpolationOrder;
66
67 TEUCHOS_TEST_FOR_EXCEPTION((LFineNodesPerDir.size() != 3) && (LFineNodesPerDir.size() != numDimensions),
69 "LFineNodesPerDir has to be of size 3 or of size numDimensions!");
70
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);
75
76 // Load coarse rate, being careful about formating
77 // Also load lFineNodesPerDir
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];
85 }
86 } else {
87 lFineNodesPerDir_h(dim) = 1;
88 coarseRate_h(dim) = 1;
89 }
90 }
91
92 Kokkos::deep_copy(lFineNodesPerDir, lFineNodesPerDir_h);
93 Kokkos::deep_copy(coarseRate, coarseRate_h);
94
95} // setupIM
96
97template <class LocalOrdinal, class GlobalOrdinal, class Node>
99 RCP<Teuchos::FancyOStream> out;
100 if (const char* dbg = std::getenv("MUELU_INDEXMANAGER_DEBUG")) {
101 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
102 out->setShowAllFrontMatter(false).setShowProcRank(true);
103 } else {
104 out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
105 }
106
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);
109
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);
114
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);
122 }
123
124 } else { // Default value for dim >= numDimensions
125 endRate_h(dim) = 1;
126 }
127 }
128
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;
133
134 // Here one element can represent either the degenerate case of one node or the more general
135 // case of two nodes, i.e. x---x is a 1D element with two nodes and x is a 1D element with
136 // one node. This helps generating a 3D space from tensorial products...
137 // A good way to handle this would be to generalize the algorithm to take into account the
138 // discretization order used in each direction, at least in the FEM sense, since a 0 degree
139 // discretization will have a unique node per element. This way 1D discretization can be
140 // viewed as a 3D problem with one 0 degree element in the y direction and one 0 degre
141 // element in the z direction.
142 // !!! Operations below are aftecting both local and global values that have two !!!
143 // different orientations. Orientations can be interchanged using mapDirG2L and mapDirL2G.
144 // coarseRate, endRate and offsets are in the global basis, as well as all the variables
145 // starting with a g.
146 // !!! while the variables starting with an l are in the local basis. !!!
147 for (int dim = 0; dim < 3; ++dim) {
148 if (dim < numDimensions) {
149 // Check whether the partition includes the "end" of the mesh which means that endRate
150 // will apply. Also make sure that endRate is not 0 which means that the mesh does not
151 // require a particular treatment at the boundaries.
152 coarseNodesPerDir_h(dim) = (lFineNodesPerDir_h(dim) - endRate_h(dim) - 1) / coarseRate_h(dim) + 2;
153
154 } else { // Default value for dim >= numDimensions
155 // endRate[dim] = 1;
156 coarseNodesPerDir_h(dim) = 1;
157 } // if (dim < numDimensions)
158
159 // This would happen if the rank does not own any nodes but in that case a subcommunicator
160 // should be used so this should really not be a concern.
161 if (lFineNodesPerDir_h(dim) < 1) {
162 coarseNodesPerDir_h(dim) = 0;
163 }
164 } // Loop for dim=0:3
165
166 // Compute cummulative values
167 numCoarseNodes10 = coarseNodesPerDir_h(0) * coarseNodesPerDir_h(1);
168 numCoarseNodes = numCoarseNodes10 * coarseNodesPerDir_h(2);
169
170 *out << "coarseNodesPerDir: {" << coarseNodesPerDir_h(0) << ", "
171 << coarseNodesPerDir_h(1) << ", " << coarseNodesPerDir_h(2) << "}" << std::endl;
172 *out << "numCoarseNodes=" << numCoarseNodes << std::endl;
173
174 // Copy Host data to Device.
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);
179}
180
181template <class LocalOrdinal, class GlobalOrdinal, class Node>
184 typename LOTupleView::HostMirror coarseNodesPerDir_h = Kokkos::create_mirror_view(coarseNodesPerDir);
185 Kokkos::deep_copy(coarseNodesPerDir_h, coarseNodesPerDir);
186 Array<LO> coarseNodesPerDirArray(3);
187
188 for (int dim = 0; dim < 3; ++dim) {
189 coarseNodesPerDirArray[dim] = coarseNodesPerDir_h(dim);
190 }
191
192 return coarseNodesPerDirArray;
193} // getCoarseNodesData
194
195} // namespace MueLu
196
197#define MUELU_INDEXMANAGER_KOKKOS_SHORT
198#endif // MUELU_INDEXMANAGER_DEF_KOKKOS_HPP
Exception throws to report errors in the internal logical of the program.
LO numCoarseNodes10
local number of nodes per 0-1 slice remaining after coarsening.
LOTupleView coarseNodesPerDir
local number of nodes per direction remaing after coarsening.
intTupleView endRate
adapted coarsening rate at the edge of the mesh in each direction.
LO lNumFineNodes10
local number of nodes per 0-1 slice.
LOTupleView lFineNodesPerDir
local number of nodes per direction.
LO numCoarseNodes
local number of nodes remaining after coarsening.
IndexManager_kokkos()=default
Default constructor, return empty object.
intTupleView coarseRate
coarsening rate in each direction
int numDimensions
Number of spacial dimensions in the problem.
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.
Namespace for MueLu classes and methods.