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// ***********************************************************************
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_DEF_KOKKOS_HPP
47#define MUELU_INDEXMANAGER_DEF_KOKKOS_HPP
48
49#include <utility>
50
51#include "Teuchos_OrdinalTraits.hpp"
52
53#include <Xpetra_MapFactory.hpp>
54
55#include "MueLu_ConfigDefs.hpp"
56#include "MueLu_Types.hpp"
57#include "MueLu_Exceptions.hpp"
59
60/*****************************************************************************
61
62****************************************************************************/
63
64namespace MueLu {
65
66 template<class LocalOrdinal, class GlobalOrdinal, class Node>
68 IndexManager_kokkos(const int NumDimensions,
69 const int interpolationOrder,
70 const int MyRank,
71 const ArrayView<const LO> LFineNodesPerDir,
72 const ArrayView<const int> CoarseRate) :
73 myRank(MyRank), coarseRate("coarsening rate"), endRate("endRate"),
74 lFineNodesPerDir("lFineNodesPerDir"), coarseNodesPerDir("lFineNodesPerDir") {
75
76 RCP<Teuchos::FancyOStream> out;
77 if(const char* dbg = std::getenv("MUELU_INDEXMANAGER_DEBUG")) {
78 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
79 out->setShowAllFrontMatter(false).setShowProcRank(true);
80 } else {
81 out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
82 }
83
84 setupIM(NumDimensions, interpolationOrder, CoarseRate, LFineNodesPerDir);
85
86 *out << "Done setting up the IndexManager" << std::endl;
87
89
90 *out << "Computed Mesh Parameters" << std::endl;
91
92 } // IndexManager_kokkos Constructor
93
94 template<class LocalOrdinal, class GlobalOrdinal, class Node>
96 setupIM(const int NumDimensions, const int interpolationOrder,
97 const ArrayView<const int> CoarseRate, const ArrayView<const LO> LFineNodesPerDir) {
98
99 numDimensions = NumDimensions;
100 interpolationOrder_ = interpolationOrder;
101
102 TEUCHOS_TEST_FOR_EXCEPTION((LFineNodesPerDir.size() != 3)
103 && (LFineNodesPerDir.size() != numDimensions),
105 "LFineNodesPerDir has to be of size 3 or of size numDimensions!");
106
107 typename Kokkos::View<LO[3], device_type>::HostMirror lFineNodesPerDir_h = Kokkos::create_mirror_view(lFineNodesPerDir);
108 Kokkos::deep_copy(lFineNodesPerDir_h, lFineNodesPerDir);
109 typename Kokkos::View<LO[3], device_type>::HostMirror coarseRate_h = Kokkos::create_mirror_view(coarseRate);
110 Kokkos::deep_copy(coarseRate_h, coarseRate);
111
112 // Load coarse rate, being careful about formating
113 // Also load lFineNodesPerDir
114 for(int dim = 0; dim < 3; ++dim) {
115 if(dim < getNumDimensions()) {
116 lFineNodesPerDir_h(dim) = LFineNodesPerDir[dim];
117 if(CoarseRate.size() == 1) {
118 coarseRate_h(dim) = CoarseRate[0];
119 } else if(CoarseRate.size() == getNumDimensions()) {
120 coarseRate_h(dim) = CoarseRate[dim];
121 }
122 } else {
123 lFineNodesPerDir_h(dim) = 1;
124 coarseRate_h(dim) = 1;
125 }
126 }
127
128 Kokkos::deep_copy(lFineNodesPerDir, lFineNodesPerDir_h);
129 Kokkos::deep_copy(coarseRate, coarseRate_h);
130
131 } // setupIM
132
133 template<class LocalOrdinal, class GlobalOrdinal, class Node>
135
136 RCP<Teuchos::FancyOStream> out;
137 if(const char* dbg = std::getenv("MUELU_INDEXMANAGER_DEBUG")) {
138 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
139 out->setShowAllFrontMatter(false).setShowProcRank(true);
140 } else {
141 out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
142 }
143
144 typename Kokkos::View<int[3], device_type>::HostMirror coarseRate_h = Kokkos::create_mirror_view(coarseRate);
145 typename Kokkos::View<int[3], device_type>::HostMirror endRate_h = Kokkos::create_mirror_view(endRate);
146
147
148 typename Kokkos::View<LO[3], device_type>::HostMirror lFineNodesPerDir_h = Kokkos::create_mirror_view(lFineNodesPerDir);
149 typename Kokkos::View<LO[3], device_type>::HostMirror coarseNodesPerDir_h = Kokkos::create_mirror_view(coarseNodesPerDir);
150 Kokkos::deep_copy(lFineNodesPerDir_h, lFineNodesPerDir);
151 Kokkos::deep_copy(coarseRate_h, coarseRate);
152
153 lNumFineNodes10 = lFineNodesPerDir_h(1)*lFineNodesPerDir_h(0);
154 lNumFineNodes = lFineNodesPerDir_h(2)*lNumFineNodes10;
155 for(int dim = 0; dim < 3; ++dim) {
156 if(dim < numDimensions) {
157 endRate_h(dim) = (lFineNodesPerDir_h(dim) - 1) % coarseRate_h(dim);
158 if(endRate_h(dim) == 0) {endRate_h(dim) = coarseRate_h(dim);}
159
160 } else { // Default value for dim >= numDimensions
161 endRate_h(dim) = 1;
162 }
163 }
164
165 *out << "lFineNodesPerDir: {" << lFineNodesPerDir_h(0) << ", " << lFineNodesPerDir_h(1) << ", "
166 << lFineNodesPerDir_h(2) << "}" << std::endl;
167 *out << "endRate: {" << endRate_h(0) << ", " << endRate_h(1) << ", "
168 << endRate_h(2) << "}" << std::endl;
169
170 // Here one element can represent either the degenerate case of one node or the more general
171 // case of two nodes, i.e. x---x is a 1D element with two nodes and x is a 1D element with
172 // one node. This helps generating a 3D space from tensorial products...
173 // A good way to handle this would be to generalize the algorithm to take into account the
174 // discretization order used in each direction, at least in the FEM sense, since a 0 degree
175 // discretization will have a unique node per element. This way 1D discretization can be
176 // viewed as a 3D problem with one 0 degree element in the y direction and one 0 degre
177 // element in the z direction.
178 // !!! Operations below are aftecting both local and global values that have two !!!
179 // different orientations. Orientations can be interchanged using mapDirG2L and mapDirL2G.
180 // coarseRate, endRate and offsets are in the global basis, as well as all the variables
181 // starting with a g.
182 // !!! while the variables starting with an l are in the local basis. !!!
183 for(int dim = 0; dim < 3; ++dim) {
184 if(dim < numDimensions) {
185 // Check whether the partition includes the "end" of the mesh which means that endRate
186 // will apply. Also make sure that endRate is not 0 which means that the mesh does not
187 // require a particular treatment at the boundaries.
188 coarseNodesPerDir_h(dim) = (lFineNodesPerDir_h(dim) - endRate_h(dim) - 1)
189 / coarseRate_h(dim) + 2;
190
191 } else { // Default value for dim >= numDimensions
192 // endRate[dim] = 1;
193 coarseNodesPerDir_h(dim) = 1;
194 } // if (dim < numDimensions)
195
196 // This would happen if the rank does not own any nodes but in that case a subcommunicator
197 // should be used so this should really not be a concern.
198 if(lFineNodesPerDir_h(dim) < 1) {coarseNodesPerDir_h(dim) = 0;}
199 } // Loop for dim=0:3
200
201 // Compute cummulative values
202 numCoarseNodes10 = coarseNodesPerDir_h(0)*coarseNodesPerDir_h(1);
203 numCoarseNodes = numCoarseNodes10*coarseNodesPerDir_h(2);
204
205 *out << "coarseNodesPerDir: {" << coarseNodesPerDir_h(0) << ", "
206 << coarseNodesPerDir_h(1) << ", " << coarseNodesPerDir_h(2) << "}" << std::endl;
207 *out << "numCoarseNodes=" << numCoarseNodes << std::endl;
208
209 // Copy Host data to Device.
210 Kokkos::deep_copy(coarseRate, coarseRate_h);
211 Kokkos::deep_copy(endRate, endRate_h);
212 Kokkos::deep_copy(lFineNodesPerDir, lFineNodesPerDir_h);
213 Kokkos::deep_copy(coarseNodesPerDir, coarseNodesPerDir_h);
214 }
215
216 template<class LocalOrdinal, class GlobalOrdinal, class Node>
219 typename LOTupleView::HostMirror coarseNodesPerDir_h = Kokkos::create_mirror_view(coarseNodesPerDir);
220 Kokkos::deep_copy(coarseNodesPerDir_h, coarseNodesPerDir);
221 Array<LO> coarseNodesPerDirArray(3);
222
223 for(int dim = 0; dim < 3; ++dim) {
224 coarseNodesPerDirArray[dim] = coarseNodesPerDir_h(dim);
225 }
226
227 return coarseNodesPerDirArray;
228 } // getCoarseNodesData
229
230} //namespace MueLu
231
232#define MUELU_INDEXMANAGER_KOKKOS_SHORT
233#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.