MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_DistanceLaplacianDropping.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_DISTANCELAPLACIANDROPPING_HPP
11#define MUELU_DISTANCELAPLACIANDROPPING_HPP
12
13#include "Kokkos_Macros.hpp"
14#include "KokkosBatched_LU_Decl.hpp"
15#include "KokkosBatched_LU_Serial_Impl.hpp"
16#include "KokkosBatched_Trsv_Decl.hpp"
17#include "KokkosBatched_Trsv_Serial_Impl.hpp"
19#include "Kokkos_Core.hpp"
20#include "Kokkos_ArithTraits.hpp"
21#include "Teuchos_RCP.hpp"
22#include "Xpetra_Matrix.hpp"
23#include "Xpetra_MultiVector.hpp"
24#include "Xpetra_MultiVectorFactory.hpp"
25
27
32template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
34 private:
35 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
36 using local_matrix_type = typename matrix_type::local_matrix_type;
37 using scalar_type = typename local_matrix_type::value_type;
39 using ATS = Kokkos::ArithTraits<scalar_type>;
40 using impl_scalar_type = typename ATS::val_type;
41 using implATS = Kokkos::ArithTraits<impl_scalar_type>;
42 using magnitudeType = typename implATS::magnitudeType;
43 using magATS = Kokkos::ArithTraits<magnitudeType>;
44 using coords_type = Xpetra::MultiVector<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>;
45 using local_coords_type = typename coords_type::dual_view_type_const::t_dev;
46
47 Teuchos::RCP<coords_type> coordsMV;
48 Teuchos::RCP<coords_type> ghostedCoordsMV;
49
52
53 public:
54 UnweightedDistanceFunctor(matrix_type& A, Teuchos::RCP<coords_type>& coords_) {
55 coordsMV = coords_;
56 auto importer = A.getCrsGraph()->getImporter();
57 if (!importer.is_null()) {
58 ghostedCoordsMV = Xpetra::MultiVectorFactory<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), coordsMV->getNumVectors());
59 ghostedCoordsMV->doImport(*coordsMV, *importer, Xpetra::INSERT);
60 coords = coordsMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
61 ghostedCoords = ghostedCoordsMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
62 } else {
63 coords = coordsMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
65 }
66 }
67
68 KOKKOS_FORCEINLINE_FUNCTION
70 magnitudeType d = magATS::zero();
72 for (size_t j = 0; j < coords.extent(1); ++j) {
73 s = coords(row, j) - ghostedCoords(col, j);
74 d += s * s;
75 }
76 return d;
77 }
78};
79
80template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
82 private:
83 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
84 using local_matrix_type = typename matrix_type::local_matrix_type;
85 using scalar_type = typename local_matrix_type::value_type;
87 using ATS = Kokkos::ArithTraits<scalar_type>;
88 using impl_scalar_type = typename ATS::val_type;
89 using implATS = Kokkos::ArithTraits<impl_scalar_type>;
90 using magnitudeType = typename implATS::magnitudeType;
91 using magATS = Kokkos::ArithTraits<magnitudeType>;
92 using coords_type = Xpetra::MultiVector<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>;
93 using local_coords_type = typename coords_type::dual_view_type_const::t_dev;
94 using material_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
95 using local_material_type = typename material_type::dual_view_type_const::t_dev;
96
97 Teuchos::RCP<coords_type> coordsMV;
98 Teuchos::RCP<coords_type> ghostedCoordsMV;
99
102
103 Teuchos::RCP<material_type> materialMV;
104 Teuchos::RCP<material_type> ghostedMaterialMV;
105
108
109 public:
110 ScalarMaterialDistanceFunctor(matrix_type& A, Teuchos::RCP<coords_type>& coords_, Teuchos::RCP<material_type>& material_) {
111 coordsMV = coords_;
112 materialMV = material_;
113 auto importer = A.getCrsGraph()->getImporter();
114 if (!importer.is_null()) {
115 ghostedCoordsMV = Xpetra::MultiVectorFactory<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), coordsMV->getNumVectors());
116 ghostedCoordsMV->doImport(*coordsMV, *importer, Xpetra::INSERT);
117 coords = coordsMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
118 ghostedCoords = ghostedCoordsMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
119
120 ghostedMaterialMV = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), materialMV->getNumVectors());
121 ghostedMaterialMV->doImport(*materialMV, *importer, Xpetra::INSERT);
122 material = materialMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
123 ghostedMaterial = ghostedMaterialMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
124 } else {
125 coords = coordsMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
127
128 material = materialMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
130 }
131 }
132
133 KOKKOS_INLINE_FUNCTION
135 // || x_row - x_col ||_S^2
136 // where
137 // S = 1/material(row) * Identity
138 magnitudeType d = magATS::zero();
139 magnitudeType d_row = magATS::zero();
140 magnitudeType d_col = magATS::zero();
142
143 for (size_t j = 0; j < coords.extent(1); ++j) {
144 s = coords(row, j) - ghostedCoords(col, j);
145 d += s * s;
146 }
147
148 d_row = d / implATS::magnitude(ghostedMaterial(row, 0));
149 d_col = d / implATS::magnitude(ghostedMaterial(col, 0));
150
151 return Kokkos::max(d_row, d_col);
152 }
153};
154
155template <class local_ordinal_type, class material_vector_type, class material_matrix_type>
157 private:
158 material_vector_type material_vector;
159 material_matrix_type material_matrix;
160
161 public:
162 TensorInversion(material_vector_type& material_vector_, material_matrix_type& material_matrix_)
163 : material_vector(material_vector_)
164 , material_matrix(material_matrix_) {}
165
166 KOKKOS_FORCEINLINE_FUNCTION
167 void operator()(local_ordinal_type i) const {
168 for (size_t j = 0; j < material_matrix.extent(1); ++j) {
169 for (size_t k = 0; k < material_matrix.extent(2); ++k) {
170 material_matrix(i, j, k) = material_vector(i, j * material_matrix.extent(1) + k);
171 }
172 }
173 auto matrix_material = Kokkos::subview(material_matrix, i, Kokkos::ALL(), Kokkos::ALL());
174 KokkosBatched::SerialLU<KokkosBatched::Algo::LU::Unblocked>::invoke(matrix_material);
175 }
176};
177
178template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
180 private:
181 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
182 using local_matrix_type = typename matrix_type::local_matrix_type;
183 using scalar_type = typename local_matrix_type::value_type;
185 using ATS = Kokkos::ArithTraits<scalar_type>;
186 using impl_scalar_type = typename ATS::val_type;
187 using implATS = Kokkos::ArithTraits<impl_scalar_type>;
188 using magnitudeType = typename implATS::magnitudeType;
189 using magATS = Kokkos::ArithTraits<magnitudeType>;
190 using coords_type = Xpetra::MultiVector<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>;
191 using local_coords_type = typename coords_type::dual_view_type_const::t_dev;
192 using material_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
193 using memory_space = typename local_matrix_type::memory_space;
194
195 using local_material_type = Kokkos::View<impl_scalar_type***, memory_space>;
196 using local_dist_type = Kokkos::View<impl_scalar_type**, memory_space>;
197
198 Teuchos::RCP<coords_type> coordsMV;
199 Teuchos::RCP<coords_type> ghostedCoordsMV;
200
203
205
207
208 const scalar_type one = ATS::one();
209
210 public:
211 TensorMaterialDistanceFunctor(matrix_type& A, Teuchos::RCP<coords_type>& coords_, Teuchos::RCP<material_type>& material_) {
212 coordsMV = coords_;
213
214 auto importer = A.getCrsGraph()->getImporter();
215 if (!importer.is_null()) {
216 ghostedCoordsMV = Xpetra::MultiVectorFactory<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), coordsMV->getNumVectors());
217 ghostedCoordsMV->doImport(*coordsMV, *importer, Xpetra::INSERT);
218 coords = coordsMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
219 ghostedCoords = ghostedCoordsMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
220 } else {
221 coords = coordsMV->getDeviceLocalView(Xpetra::Access::ReadOnly);
223 }
224
225 {
226 Teuchos::RCP<material_type> ghostedMaterial;
227 if (!importer.is_null()) {
228 ghostedMaterial = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), material_->getNumVectors());
229 ghostedMaterial->doImport(*material_, *importer, Xpetra::INSERT);
230 } else {
231 ghostedMaterial = material_;
232 }
233
234 using execution_space = typename Node::execution_space;
235 using range_type = Kokkos::RangePolicy<LocalOrdinal, execution_space>;
236
237 local_ordinal_type dim = std::sqrt(material_->getNumVectors());
238 auto lclMaterial = ghostedMaterial->getDeviceLocalView(Xpetra::Access::ReadOnly);
239 material = local_material_type("material", lclMaterial.extent(0), dim, dim);
240 lcl_dist = local_dist_type("material", lclMaterial.extent(0), dim);
242 Kokkos::parallel_for("MueLu:TensorMaterialDistanceFunctor::inversion", range_type(0, lclMaterial.extent(0)), functor);
243 }
244 }
245
246 KOKKOS_INLINE_FUNCTION
248 // || x_row - x_col ||_S^2
249 // where
250 // S = inv(material(col))
251
252 // row material
253 impl_scalar_type d_row = implATS::zero();
254 {
255 auto matrix_row_material = Kokkos::subview(material, row, Kokkos::ALL(), Kokkos::ALL());
256 auto dist = Kokkos::subview(lcl_dist, row, Kokkos::ALL());
257
258 for (size_t j = 0; j < coords.extent(1); ++j) {
259 dist(j) = coords(row, j) - ghostedCoords(col, j);
260 }
261
262 KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Lower, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::Unit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(one, matrix_row_material, dist);
263 KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Upper, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::NonUnit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(one, matrix_row_material, dist);
264
265 for (size_t j = 0; j < coords.extent(1); ++j) {
266 d_row += dist(j) * (coords(row, j) - ghostedCoords(col, j));
267 }
268 }
269
270 // column material
271 impl_scalar_type d_col = implATS::zero();
272 {
273 auto matrix_col_material = Kokkos::subview(material, col, Kokkos::ALL(), Kokkos::ALL());
274 auto dist = Kokkos::subview(lcl_dist, row, Kokkos::ALL());
275
276 for (size_t j = 0; j < coords.extent(1); ++j) {
277 dist(j) = coords(row, j) - ghostedCoords(col, j);
278 }
279
280 KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Lower, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::Unit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(one, matrix_col_material, dist);
281 KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Upper, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::NonUnit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(one, matrix_col_material, dist);
282
283 for (size_t j = 0; j < coords.extent(1); ++j) {
284 d_col += dist(j) * (coords(row, j) - ghostedCoords(col, j));
285 }
286 }
287
288 return Kokkos::max(implATS::magnitude(d_row), implATS::magnitude(d_col));
289 }
290};
291
295template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class DistanceFunctorType>
296Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
297getDiagonal(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
298 DistanceFunctorType& distFunctor) {
299 using scalar_type = Scalar;
300 using local_ordinal_type = LocalOrdinal;
301 using global_ordinal_type = GlobalOrdinal;
302 using node_type = Node;
303 using ATS = Kokkos::ArithTraits<scalar_type>;
304 using impl_scalar_type = typename ATS::val_type;
305 using implATS = Kokkos::ArithTraits<impl_scalar_type>;
306 using magnitudeType = typename implATS::magnitudeType;
307 using execution_space = typename Node::execution_space;
308 using range_type = Kokkos::RangePolicy<LocalOrdinal, execution_space>;
309
310 auto diag = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A.getRowMap(), 1);
311 {
312 auto lclA = A.getLocalMatrixDevice();
313 auto lclDiag = diag->getDeviceLocalView(Xpetra::Access::OverwriteAll);
314
315 Kokkos::parallel_for(
316 "MueLu:CoalesceDropF:Build:scalar_filter:laplacian_diag",
317 range_type(0, lclA.numRows()),
318 KOKKOS_LAMBDA(const local_ordinal_type& row) {
319 auto rowView = lclA.rowConst(row);
320 auto length = rowView.length;
321
322 magnitudeType d;
323 impl_scalar_type d2 = implATS::zero();
324 for (local_ordinal_type colID = 0; colID < length; colID++) {
325 auto col = rowView.colidx(colID);
326 if (row != col) {
327 d = distFunctor.distance2(row, col);
328 d2 += implATS::one() / d;
329 }
330 }
331 lclDiag(row, 0) = d2;
332 });
333 }
334 auto importer = A.getCrsGraph()->getImporter();
335 if (!importer.is_null()) {
336 auto ghostedDiag = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A.getColMap(), 1);
337 ghostedDiag->doImport(*diag, *importer, Xpetra::INSERT);
338 return ghostedDiag;
339 } else {
340 return diag;
341 }
342}
343
354template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class DistanceFunctorType>
356 private:
357 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
358 using local_matrix_type = typename matrix_type::local_matrix_type;
359 using scalar_type = typename local_matrix_type::value_type;
360 using local_ordinal_type = typename local_matrix_type::ordinal_type;
361 using memory_space = typename local_matrix_type::memory_space;
362 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
363 using diag_view_type = typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
364
365 using results_view = Kokkos::View<DecisionType*, memory_space>;
366
367 using ATS = Kokkos::ArithTraits<scalar_type>;
368 using magnitudeType = typename ATS::magnitudeType;
369 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
370
373 Teuchos::RCP<diag_vec_type> diagVec;
374 diag_view_type diag; // corresponds to overlapped diagonal
375 DistanceFunctorType dist2;
377 const scalar_type one = ATS::one();
378
379 public:
380 DropFunctor(matrix_type& A_, magnitudeType threshold, DistanceFunctorType& dist2_, results_view& results_)
381 : A(A_.getLocalMatrixDevice())
382 , eps(threshold)
383 , dist2(dist2_)
384 , results(results_) {
386 auto lclDiag2d = diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
387 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
388 }
389
390 KOKKOS_FORCEINLINE_FUNCTION
392 auto row = A.rowConst(rlid);
393 const size_t offset = A.graph.row_map(rlid);
394 for (local_ordinal_type k = 0; k < row.length; ++k) {
395 auto clid = row.colidx(k);
396
397 scalar_type val;
398 if (rlid != clid) {
399 val = one / dist2.distance2(rlid, clid);
400 } else {
401 val = diag(rlid);
402 }
403 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
404 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
405
406 results(offset + k) = Kokkos::max((aij2 <= eps * eps * aiiajj) ? DROP : KEEP,
407 results(offset + k));
408 }
409 }
410};
411
412} // namespace MueLu::DistanceLaplacian
413
414#endif
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
typename local_matrix_type::value_type scalar_type
Kokkos::View< DecisionType *, memory_space > results_view
typename Kokkos::DualView< const scalar_type *, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged >::t_dev diag_view_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
typename matrix_type::local_matrix_type local_matrix_type
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< const bool *, memory_space > boundary_nodes_view
typename local_matrix_type::memory_space memory_space
DropFunctor(matrix_type &A_, magnitudeType threshold, DistanceFunctorType &dist2_, results_view &results_)
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > material_type
Xpetra::MultiVector< magnitudeType, LocalOrdinal, GlobalOrdinal, Node > coords_type
KOKKOS_INLINE_FUNCTION magnitudeType distance2(const local_ordinal_type row, const local_ordinal_type col) const
typename coords_type::dual_view_type_const::t_dev local_coords_type
ScalarMaterialDistanceFunctor(matrix_type &A, Teuchos::RCP< coords_type > &coords_, Teuchos::RCP< material_type > &material_)
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
typename material_type::dual_view_type_const::t_dev local_material_type
TensorInversion(material_vector_type &material_vector_, material_matrix_type &material_matrix_)
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type i) const
Kokkos::View< impl_scalar_type **, memory_space > local_dist_type
typename coords_type::dual_view_type_const::t_dev local_coords_type
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > material_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
Kokkos::View< impl_scalar_type ***, memory_space > local_material_type
Xpetra::MultiVector< magnitudeType, LocalOrdinal, GlobalOrdinal, Node > coords_type
KOKKOS_INLINE_FUNCTION magnitudeType distance2(const local_ordinal_type row, const local_ordinal_type col) const
TensorMaterialDistanceFunctor(matrix_type &A, Teuchos::RCP< coords_type > &coords_, Teuchos::RCP< material_type > &material_)
KOKKOS_FORCEINLINE_FUNCTION magnitudeType distance2(const local_ordinal_type row, const local_ordinal_type col) const
Xpetra::MultiVector< magnitudeType, LocalOrdinal, GlobalOrdinal, Node > coords_type
typename coords_type::dual_view_type_const::t_dev local_coords_type
UnweightedDistanceFunctor(matrix_type &A, Teuchos::RCP< coords_type > &coords_)
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getDiagonal(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, DistanceFunctorType &distFunctor)