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 "KokkosKernels_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_device_type;
37 using scalar_type = typename local_matrix_type::value_type;
39 using ATS = KokkosKernels::ArithTraits<scalar_type>;
40 using impl_scalar_type = typename ATS::val_type;
41 using implATS = KokkosKernels::ArithTraits<impl_scalar_type>;
42 using magnitudeType = typename implATS::magnitudeType;
43 using magATS = KokkosKernels::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(), false);
59 ghostedCoordsMV->doImport(*coordsMV, *importer, Xpetra::INSERT);
60 coords = coordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
61 ghostedCoords = ghostedCoordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
62 } else {
63 coords = coordsMV->getLocalViewDevice(Tpetra::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
84template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class weight_type>
86 private:
87 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
88 using local_matrix_type = typename matrix_type::local_matrix_device_type;
89 using scalar_type = typename local_matrix_type::value_type;
91 using ATS = KokkosKernels::ArithTraits<scalar_type>;
92 using impl_scalar_type = typename ATS::val_type;
93 using implATS = KokkosKernels::ArithTraits<impl_scalar_type>;
94 using magnitudeType = typename implATS::magnitudeType;
95 using magATS = KokkosKernels::ArithTraits<magnitudeType>;
96 using coords_type = Xpetra::MultiVector<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>;
97 using local_coords_type = typename coords_type::dual_view_type_const::t_dev;
98
99 Teuchos::RCP<coords_type> coordsMV;
100 Teuchos::RCP<coords_type> ghostedCoordsMV;
101
104
105 weight_type weight;
106
107 public:
108 WeightedDistanceFunctor(matrix_type& A, Teuchos::RCP<coords_type>& coords_, weight_type weight_) {
109 coordsMV = coords_;
110 auto importer = A.getCrsGraph()->getImporter();
111 if (!importer.is_null()) {
112 ghostedCoordsMV = Xpetra::MultiVectorFactory<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), coordsMV->getNumVectors(), false);
113 ghostedCoordsMV->doImport(*coordsMV, *importer, Xpetra::INSERT);
114 coords = coordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
115 ghostedCoords = ghostedCoordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
116 } else {
117 coords = coordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
119 }
120 weight = weight_;
121 TEUCHOS_ASSERT(weight.extent(0) == coordsMV->getNumVectors());
122 }
123
124 KOKKOS_FORCEINLINE_FUNCTION
126 magnitudeType d = magATS::zero();
129 for (size_t j = 0; j < coords.extent(1); ++j) {
130 s = coords(row, j) - ghostedCoords(col, j);
131 w = weight(j);
132 d += w * s * s;
133 }
134 return d;
135 }
136};
137
142template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class weight_type>
144 private:
145 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
146 using local_matrix_type = typename matrix_type::local_matrix_device_type;
147 using scalar_type = typename local_matrix_type::value_type;
149 using ATS = KokkosKernels::ArithTraits<scalar_type>;
150 using impl_scalar_type = typename ATS::val_type;
151 using implATS = KokkosKernels::ArithTraits<impl_scalar_type>;
152 using magnitudeType = typename implATS::magnitudeType;
153 using magATS = KokkosKernels::ArithTraits<magnitudeType>;
154 using coords_type = Xpetra::MultiVector<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>;
155 using local_coords_type = typename coords_type::dual_view_type_const::t_dev;
156
157 Teuchos::RCP<coords_type> coordsMV;
158 Teuchos::RCP<coords_type> ghostedCoordsMV;
159
162
163 weight_type weight;
165
166 public:
167 BlockWeightedDistanceFunctor(matrix_type& A, Teuchos::RCP<coords_type>& coords_, weight_type weight_, local_ordinal_type interleaved_blocksize_) {
168 coordsMV = coords_;
169 auto importer = A.getCrsGraph()->getImporter();
170 if (!importer.is_null()) {
171 ghostedCoordsMV = Xpetra::MultiVectorFactory<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), coordsMV->getNumVectors(), false);
172 ghostedCoordsMV->doImport(*coordsMV, *importer, Xpetra::INSERT);
173 coords = coordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
174 ghostedCoords = ghostedCoordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
175 } else {
176 coords = coordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
178 }
179 weight = weight_;
180 interleaved_blocksize = interleaved_blocksize_;
181 TEUCHOS_ASSERT(weight.extent(0) == coordsMV->getNumVectors());
182 }
183
184 KOKKOS_FORCEINLINE_FUNCTION
186 magnitudeType d = magATS::zero();
190 local_ordinal_type block_start = block_id * interleaved_blocksize;
191 for (size_t j = 0; j < coords.extent(1); ++j) {
192 s = coords(row, j) - ghostedCoords(col, j);
193 w = weight(block_start + j);
194 d += w * s * s;
195 }
196 return d;
197 }
198};
199
200template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
202 private:
203 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
204 using local_matrix_type = typename matrix_type::local_matrix_device_type;
205 using scalar_type = typename local_matrix_type::value_type;
207 using ATS = KokkosKernels::ArithTraits<scalar_type>;
208 using impl_scalar_type = typename ATS::val_type;
209 using implATS = KokkosKernels::ArithTraits<impl_scalar_type>;
210 using magnitudeType = typename implATS::magnitudeType;
211 using magATS = KokkosKernels::ArithTraits<magnitudeType>;
212 using coords_type = Xpetra::MultiVector<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>;
213 using local_coords_type = typename coords_type::dual_view_type_const::t_dev;
214 using material_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
215 using local_material_type = typename material_type::dual_view_type_const::t_dev;
216
217 Teuchos::RCP<coords_type> coordsMV;
218 Teuchos::RCP<coords_type> ghostedCoordsMV;
219
222
223 Teuchos::RCP<material_type> materialMV;
224 Teuchos::RCP<material_type> ghostedMaterialMV;
225
228
229 public:
230 ScalarMaterialDistanceFunctor(matrix_type& A, Teuchos::RCP<coords_type>& coords_, Teuchos::RCP<material_type>& material_) {
231 coordsMV = coords_;
232 materialMV = material_;
233 auto importer = A.getCrsGraph()->getImporter();
234 if (!importer.is_null()) {
235 ghostedCoordsMV = Xpetra::MultiVectorFactory<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), coordsMV->getNumVectors(), false);
236 ghostedCoordsMV->doImport(*coordsMV, *importer, Xpetra::INSERT);
237 coords = coordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
238 ghostedCoords = ghostedCoordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
239
240 ghostedMaterialMV = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), materialMV->getNumVectors(), false);
241 ghostedMaterialMV->doImport(*materialMV, *importer, Xpetra::INSERT);
242 material = materialMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
243 ghostedMaterial = ghostedMaterialMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
244 } else {
245 coords = coordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
247
248 material = materialMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
250 }
251 }
252
253 KOKKOS_INLINE_FUNCTION
255 // || x_row - x_col ||_S^2
256 // where
257 // S = 1/material(row) * Identity
258 magnitudeType d = magATS::zero();
259 magnitudeType d_row = magATS::zero();
260 magnitudeType d_col = magATS::zero();
262
263 for (size_t j = 0; j < coords.extent(1); ++j) {
264 s = coords(row, j) - ghostedCoords(col, j);
265 d += s * s;
266 }
267
268 d_row = d / implATS::magnitude(ghostedMaterial(row, 0));
269 d_col = d / implATS::magnitude(ghostedMaterial(col, 0));
270
271 return Kokkos::max(d_row, d_col);
272 }
273};
274
275template <class local_ordinal_type, class material_vector_type, class material_matrix_type>
277 private:
278 material_vector_type material_vector;
279 material_matrix_type material_matrix;
280
281 public:
282 TensorInversion(material_vector_type& material_vector_, material_matrix_type& material_matrix_)
283 : material_vector(material_vector_)
284 , material_matrix(material_matrix_) {}
285
286 KOKKOS_FORCEINLINE_FUNCTION
287 void operator()(local_ordinal_type i) const {
288 for (size_t j = 0; j < material_matrix.extent(1); ++j) {
289 for (size_t k = 0; k < material_matrix.extent(2); ++k) {
290 material_matrix(i, j, k) = material_vector(i, j * material_matrix.extent(1) + k);
291 }
292 }
293 auto matrix_material = Kokkos::subview(material_matrix, i, Kokkos::ALL(), Kokkos::ALL());
294 KokkosBatched::SerialLU<KokkosBatched::Algo::LU::Unblocked>::invoke(matrix_material);
295 }
296};
297
298template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
300 private:
301 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
302 using local_matrix_type = typename matrix_type::local_matrix_device_type;
303 using scalar_type = typename local_matrix_type::value_type;
305 using ATS = KokkosKernels::ArithTraits<scalar_type>;
306 using impl_scalar_type = typename ATS::val_type;
307 using implATS = KokkosKernels::ArithTraits<impl_scalar_type>;
308 using magnitudeType = typename implATS::magnitudeType;
309 using magATS = KokkosKernels::ArithTraits<magnitudeType>;
310 using coords_type = Xpetra::MultiVector<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>;
311 using local_coords_type = typename coords_type::dual_view_type_const::t_dev;
312 using material_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
313 using memory_space = typename local_matrix_type::memory_space;
314
315 using local_material_type = Kokkos::View<impl_scalar_type***, memory_space>;
316 using local_dist_type = Kokkos::View<impl_scalar_type**, memory_space>;
317
318 Teuchos::RCP<coords_type> coordsMV;
319 Teuchos::RCP<coords_type> ghostedCoordsMV;
320
323
325
327
328 const scalar_type one = ATS::one();
329
330 public:
331 TensorMaterialDistanceFunctor(matrix_type& A, Teuchos::RCP<coords_type>& coords_, Teuchos::RCP<material_type>& material_) {
332 coordsMV = coords_;
333
334 auto importer = A.getCrsGraph()->getImporter();
335 if (!importer.is_null()) {
336 ghostedCoordsMV = Xpetra::MultiVectorFactory<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), coordsMV->getNumVectors(), false);
337 ghostedCoordsMV->doImport(*coordsMV, *importer, Xpetra::INSERT);
338 coords = coordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
339 ghostedCoords = ghostedCoordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
340 } else {
341 coords = coordsMV->getLocalViewDevice(Tpetra::Access::ReadOnly);
343 }
344
345 {
346 Teuchos::RCP<material_type> ghostedMaterial;
347 if (!importer.is_null()) {
348 ghostedMaterial = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(), material_->getNumVectors(), false);
349 ghostedMaterial->doImport(*material_, *importer, Xpetra::INSERT);
350 } else {
351 ghostedMaterial = material_;
352 }
353
354 using execution_space = typename Node::execution_space;
355 using range_type = Kokkos::RangePolicy<LocalOrdinal, execution_space>;
356
357 local_ordinal_type dim = std::sqrt(material_->getNumVectors());
358 auto lclMaterial = ghostedMaterial->getLocalViewDevice(Tpetra::Access::ReadOnly);
359 material = local_material_type("material", lclMaterial.extent(0), dim, dim);
360 lcl_dist = local_dist_type("material", lclMaterial.extent(0), dim);
362 Kokkos::parallel_for("MueLu:TensorMaterialDistanceFunctor::inversion", range_type(0, lclMaterial.extent(0)), functor);
363 }
364 }
365
366 KOKKOS_INLINE_FUNCTION
368 // || x_row - x_col ||_S^2
369 // where
370 // S = inv(material(col))
371
372 // row material
373 impl_scalar_type d_row = implATS::zero();
374 {
375 auto matrix_row_material = Kokkos::subview(material, row, Kokkos::ALL(), Kokkos::ALL());
376 auto dist = Kokkos::subview(lcl_dist, row, Kokkos::ALL());
377
378 for (size_t j = 0; j < coords.extent(1); ++j) {
379 dist(j) = coords(row, j) - ghostedCoords(col, j);
380 }
381
382 KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Lower, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::Unit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(one, matrix_row_material, dist);
383 KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Upper, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::NonUnit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(one, matrix_row_material, dist);
384
385 for (size_t j = 0; j < coords.extent(1); ++j) {
386 d_row += dist(j) * (coords(row, j) - ghostedCoords(col, j));
387 }
388 }
389
390 // column material
391 impl_scalar_type d_col = implATS::zero();
392 {
393 auto matrix_col_material = Kokkos::subview(material, col, Kokkos::ALL(), Kokkos::ALL());
394 auto dist = Kokkos::subview(lcl_dist, row, Kokkos::ALL());
395
396 for (size_t j = 0; j < coords.extent(1); ++j) {
397 dist(j) = coords(row, j) - ghostedCoords(col, j);
398 }
399
400 KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Lower, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::Unit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(one, matrix_col_material, dist);
401 KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Upper, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::NonUnit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(one, matrix_col_material, dist);
402
403 for (size_t j = 0; j < coords.extent(1); ++j) {
404 d_col += dist(j) * (coords(row, j) - ghostedCoords(col, j));
405 }
406 }
407
408 return Kokkos::max(implATS::magnitude(d_row), implATS::magnitude(d_col));
409 }
410};
411
415template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class DistanceFunctorType>
416Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
417getDiagonal(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
418 DistanceFunctorType& distFunctor) {
419 using scalar_type = Scalar;
420 using local_ordinal_type = LocalOrdinal;
421 using ATS = KokkosKernels::ArithTraits<scalar_type>;
422 using impl_scalar_type = typename ATS::val_type;
423 using implATS = KokkosKernels::ArithTraits<impl_scalar_type>;
424 using magnitudeType = typename implATS::magnitudeType;
425 using execution_space = typename Node::execution_space;
426 using range_type = Kokkos::RangePolicy<LocalOrdinal, execution_space>;
427
428 auto diag = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A.getRowMap(), 1, false);
429 {
430 auto lclA = A.getLocalMatrixDevice();
431 auto lclDiag = diag->getLocalViewDevice(Tpetra::Access::OverwriteAll);
432
433 Kokkos::parallel_for(
434 "MueLu:CoalesceDropF:Build:scalar_filter:laplacian_diag",
435 range_type(0, lclA.numRows()),
436 KOKKOS_LAMBDA(const local_ordinal_type& row) {
437 auto rowView = lclA.rowConst(row);
438 auto length = rowView.length;
439
440 magnitudeType d;
441 impl_scalar_type d2 = implATS::zero();
442 bool haveAddedToDiag = false;
443 for (local_ordinal_type colID = 0; colID < length; colID++) {
444 auto col = rowView.colidx(colID);
445 if (row != col) {
446 d = distFunctor.distance2(row, col);
447 d2 += implATS::one() / d;
448 haveAddedToDiag = true;
449 }
450 }
451
452 // Deal with the situation where boundary conditions have only been enforced on rows, but not on columns.
453 // We enforce dropping of these entries by assigning a very large number to the diagonal entries corresponding to BCs.
454 lclDiag(row, 0) = !haveAddedToDiag ? implATS::squareroot(implATS::rmax()) : d2;
455 });
456 }
457 auto importer = A.getCrsGraph()->getImporter();
458 if (!importer.is_null()) {
459 auto ghostedDiag = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A.getColMap(), 1, false);
460 ghostedDiag->doImport(*diag, *importer, Xpetra::INSERT);
461 return ghostedDiag;
462 } else {
463 return diag;
464 }
465}
466
467template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class DistanceFunctorType>
468Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
469getMaxMinusOffDiagonal(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
470 DistanceFunctorType& distFunctor) {
471 using scalar_type = Scalar;
472 using local_ordinal_type = LocalOrdinal;
473 using ATS = KokkosKernels::ArithTraits<scalar_type>;
474 using impl_scalar_type = typename ATS::val_type;
475 using implATS = KokkosKernels::ArithTraits<impl_scalar_type>;
476 using magnitudeType = typename implATS::magnitudeType;
477 using execution_space = typename Node::execution_space;
478 using range_type = Kokkos::RangePolicy<LocalOrdinal, execution_space>;
479
480 auto diag = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A.getRowMap(), 1, false);
481 {
482 auto lclA = A.getLocalMatrixDevice();
483 auto lclDiag = diag->getLocalViewDevice(Tpetra::Access::OverwriteAll);
484
485 Kokkos::parallel_for(
486 "MueLu:CoalesceDropF:Build:scalar_filter:laplacian_diag",
487 range_type(0, lclA.numRows()),
488 KOKKOS_LAMBDA(const local_ordinal_type& row) {
489 auto rowView = lclA.rowConst(row);
490 auto length = rowView.length;
491
492 impl_scalar_type mymax = implATS::zero();
493 magnitudeType d;
494 impl_scalar_type d2;
495 for (local_ordinal_type colID = 0; colID < length; colID++) {
496 auto col = rowView.colidx(colID);
497 if (row != col) {
498 d = distFunctor.distance2(row, col);
499 d2 = implATS::one() / d;
500 if (implATS::magnitude(mymax) < implATS::magnitude(d2))
501 mymax = implATS::magnitude(d2);
502 }
503 }
504 lclDiag(row, 0) = mymax;
505 });
506 }
507 auto importer = A.getCrsGraph()->getImporter();
508 if (!importer.is_null()) {
509 auto ghostedDiag = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A.getColMap(), 1, false);
510 ghostedDiag->doImport(*diag, *importer, Xpetra::INSERT);
511 return ghostedDiag;
512 } else {
513 return diag;
514 }
515}
516
527template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class DistanceFunctorType, Misc::StrengthMeasure measure>
529 public:
530 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
531 using local_matrix_type = typename matrix_type::local_matrix_device_type;
532 using scalar_type = typename local_matrix_type::value_type;
533 using local_ordinal_type = typename local_matrix_type::ordinal_type;
534 using memory_space = typename local_matrix_type::memory_space;
535 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
536 using diag_view_type = typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
537
538 using results_view = Kokkos::View<DecisionType*, memory_space>;
539
540 using ATS = KokkosKernels::ArithTraits<scalar_type>;
541 using magnitudeType = typename ATS::magnitudeType;
542 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
543 using mATS = KokkosKernels::ArithTraits<magnitudeType>;
544
545 private:
548 Teuchos::RCP<diag_vec_type> diagVec;
549 diag_view_type diag; // corresponds to overlapped diagonal
550 DistanceFunctorType dist2;
552 const scalar_type one = ATS::one();
553
554 public:
555 DropFunctor(matrix_type& A_, magnitudeType threshold, DistanceFunctorType& dist2_, results_view& results_)
556 : A(A_.getLocalMatrixDevice())
557 , eps(threshold)
558 , dist2(dist2_)
559 , results(results_) {
560 if constexpr ((measure == Misc::SmoothedAggregationMeasure) || (measure == Misc::SignedSmoothedAggregationMeasure)) {
562 auto lclDiag2d = diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
563 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
564 } else if constexpr (measure == Misc::SignedRugeStuebenMeasure) {
566 auto lclDiag2d = diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
567 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
568 }
569 }
570
571 KOKKOS_FORCEINLINE_FUNCTION
573 auto row = A.rowConst(rlid);
574 const size_t offset = A.graph.row_map(rlid);
575
576#ifdef MUELU_COALESCE_DROP_DEBUG
577 {
578 Kokkos::printf("SoC: ");
579 for (local_ordinal_type k = 0; k < row.length; ++k) {
580 auto clid = row.colidx(k);
581
582 scalar_type val;
583 if (rlid != clid) {
584 val = -one / dist2.distance2(rlid, clid);
585 } else {
586 val = diag(rlid);
587 }
588
589 if constexpr (measure == Misc::SmoothedAggregationMeasure) {
590 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
591 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
592
593 Kokkos::printf("%5f ", ATS::sqrt(aij2 / aiiajj));
594 } else if constexpr (measure == Misc::SignedRugeStuebenMeasure) {
595 auto neg_aij = -ATS::real(val);
596 auto max_neg_aik = ATS::real(diag(rlid));
597 Kokkos::printf("%5f ", neg_aij / max_neg_aik);
598 } else if constexpr (measure == Misc::SignedSmoothedAggregationMeasure) {
599 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
600 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
601 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
602 // + |a_ij|^2, if a_ij < 0, - |a_ij|^2 if a_ij >=0
603 if (!is_nonpositive)
604 aij2 = -aij2;
605 Kokkos::printf("%5f ", ATS::sqrt(aij2 / aiiajj));
606 }
607 }
608 Kokkos::printf("\n");
609 }
610#endif
611
612 for (local_ordinal_type k = 0; k < row.length; ++k) {
613 auto clid = row.colidx(k);
614
615 scalar_type val;
616 if (rlid != clid) {
617 val = -one / dist2.distance2(rlid, clid);
618 } else {
619 val = diag(rlid);
620 }
621
622 if constexpr (measure == Misc::SmoothedAggregationMeasure) {
623 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
624 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
625
626 results(offset + k) = Kokkos::max((aij2 <= eps * eps * aiiajj) ? DROP : KEEP,
627 results(offset + k));
628 } else if constexpr (measure == Misc::SignedRugeStuebenMeasure) {
629 auto neg_aij = -ATS::real(val);
630 auto max_neg_aik = eps * ATS::real(diag(rlid));
631 results(offset + k) = Kokkos::max((neg_aij < max_neg_aik) ? DROP : KEEP,
632 results(offset + k));
633 } else if constexpr (measure == Misc::SignedSmoothedAggregationMeasure) {
634 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
635 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
636 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
637 // + |a_ij|^2, if a_ij < 0, - |a_ij|^2 if a_ij >=0
638 if (!is_nonpositive)
639 aij2 = -aij2;
640 results(offset + k) = Kokkos::max((aij2 <= eps * eps * aiiajj) ? DROP : KEEP,
641 results(offset + k));
642 }
643 }
644 }
645};
646
647template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class DistanceFunctorType, Misc::StrengthMeasure measure>
649 public:
650 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
651 using local_matrix_type = typename matrix_type::local_matrix_device_type;
652 using scalar_type = typename local_matrix_type::value_type;
653 using local_ordinal_type = typename local_matrix_type::ordinal_type;
654 using memory_space = typename local_matrix_type::memory_space;
655 using block_indices_view_type = Kokkos::View<local_ordinal_type*, memory_space>;
656 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
657 using diag_view_type = typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
658
659 using results_view = Kokkos::View<DecisionType*, memory_space>;
660
661 using ATS = KokkosKernels::ArithTraits<scalar_type>;
662 using magnitudeType = typename ATS::magnitudeType;
663 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
664 using mATS = KokkosKernels::ArithTraits<magnitudeType>;
665
666 private:
669 Teuchos::RCP<diag_vec_type> diagVec;
670 diag_view_type diag; // corresponds to overlapped diagonal
671 DistanceFunctorType dist2;
675 const scalar_type one = ATS::one();
676
677 public:
678 VectorDropFunctor(matrix_type& A_, matrix_type& mergedA_, magnitudeType threshold, DistanceFunctorType& dist2_, results_view& results_, block_indices_view_type point_to_block_, block_indices_view_type ghosted_point_to_block_)
679 : A(A_.getLocalMatrixDevice())
680 , eps(threshold)
681 , dist2(dist2_)
682 , results(results_)
683 , point_to_block(point_to_block_)
684 , ghosted_point_to_block(ghosted_point_to_block_) {
685 if constexpr ((measure == Misc::SmoothedAggregationMeasure) || (measure == Misc::SignedSmoothedAggregationMeasure)) {
686 diagVec = getDiagonal(mergedA_, dist2);
687 auto lclDiag2d = diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
688 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
689 } else if (measure == Misc::SignedRugeStuebenMeasure) {
691 auto lclDiag2d = diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
692 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
693 }
694 }
695
696 KOKKOS_FORCEINLINE_FUNCTION
698 auto brlid = point_to_block(rlid);
699 auto row = A.rowConst(rlid);
700 const size_t offset = A.graph.row_map(rlid);
701
702#ifdef MUELU_COALESCE_DROP_DEBUG
703 {
704 Kokkos::printf("SoC: ");
705 for (local_ordinal_type k = 0; k < row.length; ++k) {
706 auto clid = row.colidx(k);
707 auto bclid = ghosted_point_to_block(clid);
708
709 scalar_type val;
710 if (brlid != bclid) {
711 val = -one / dist2.distance2(brlid, bclid);
712 } else {
713 val = diag(brlid);
714 }
715
716 if constexpr (measure == Misc::SmoothedAggregationMeasure) {
717 auto aiiajj = ATS::magnitude(diag(brlid)) * ATS::magnitude(diag(bclid)); // |a_ii|*|a_jj|
718 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
719
720 Kokkos::printf("%5f ", ATS::sqrt(aij2 / aiiajj));
721 } else if constexpr (measure == Misc::SignedRugeStuebenMeasure) {
722 auto neg_aij = -ATS::real(val);
723 auto max_neg_aik = eps * ATS::real(diag(brlid));
724 Kokkos::printf("%5f ", neg_aij / max_neg_aik);
725 } else if constexpr (measure == Misc::SignedSmoothedAggregationMeasure) {
726 auto aiiajj = ATS::magnitude(diag(brlid)) * ATS::magnitude(diag(bclid)); // |a_ii|*|a_jj|
727 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
728 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
729 // + |a_ij|^2, if a_ij < 0, - |a_ij|^2 if a_ij >=0
730 if (!is_nonpositive)
731 aij2 = -aij2;
732 Kokkos::printf("%5f ", ATS::sqrt(aij2 / aiiajj));
733 }
734 }
735 Kokkos::printf("\n");
736 }
737#endif
738
739 for (local_ordinal_type k = 0; k < row.length; ++k) {
740 auto clid = row.colidx(k);
741 auto bclid = ghosted_point_to_block(clid);
742
743 scalar_type val;
744 if (brlid != bclid) {
745 val = -one / dist2.distance2(brlid, bclid);
746 } else {
747 val = diag(brlid);
748 }
749
750 if constexpr (measure == Misc::SmoothedAggregationMeasure) {
751 auto aiiajj = ATS::magnitude(diag(brlid)) * ATS::magnitude(diag(bclid)); // |a_ii|*|a_jj|
752 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
753
754 results(offset + k) = Kokkos::max((aij2 <= eps * eps * aiiajj) ? DROP : KEEP,
755 results(offset + k));
756 } else if constexpr (measure == Misc::SignedRugeStuebenMeasure) {
757 auto neg_aij = -ATS::real(val);
758 auto max_neg_aik = eps * ATS::real(diag(brlid));
759 results(offset + k) = Kokkos::max((neg_aij < max_neg_aik) ? DROP : KEEP,
760 results(offset + k));
761 } else if constexpr (measure == Misc::SignedSmoothedAggregationMeasure) {
762 auto aiiajj = ATS::magnitude(diag(brlid)) * ATS::magnitude(diag(bclid)); // |a_ii|*|a_jj|
763 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
764 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
765 // + |a_ij|^2, if a_ij < 0, - |a_ij|^2 if a_ij >=0
766 if (!is_nonpositive)
767 aij2 = -aij2;
768 results(offset + k) = Kokkos::max((aij2 <= eps * eps * aiiajj) ? DROP : KEEP,
769 results(offset + k));
770 }
771 }
772 }
773};
774
775// helper function to allow partial template deduction
776template <Misc::StrengthMeasure measure, class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class DistanceFunctorType>
777auto make_drop_functor(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A_,
779 DistanceFunctorType& dist2_,
782 return functor;
783}
784
785template <Misc::StrengthMeasure measure, class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class DistanceFunctorType>
786auto make_vector_drop_functor(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A_,
787 Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& mergedA_,
789 DistanceFunctorType& dist2_,
793 auto functor = VectorDropFunctor<Scalar, LocalOrdinal, GlobalOrdinal, Node, DistanceFunctorType, measure>(A_, mergedA_, threshold, dist2_, results_, point_to_block_, ghosted_point_to_block_);
794 return functor;
795}
796
797} // namespace MueLu::DistanceLaplacian
798
799#endif
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
BlockWeightedDistanceFunctor(matrix_type &A, Teuchos::RCP< coords_type > &coords_, weight_type weight_, local_ordinal_type interleaved_blocksize_)
typename matrix_type::local_matrix_device_type local_matrix_type
Xpetra::MultiVector< magnitudeType, LocalOrdinal, GlobalOrdinal, Node > coords_type
typename coords_type::dual_view_type_const::t_dev local_coords_type
KOKKOS_FORCEINLINE_FUNCTION magnitudeType distance2(const local_ordinal_type row, const local_ordinal_type col) const
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
Drops entries the unscaled distance Laplacian.
KokkosKernels::ArithTraits< magnitudeType > mATS
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
typename local_matrix_type::memory_space memory_space
KokkosKernels::ArithTraits< scalar_type > ATS
Kokkos::View< DecisionType *, memory_space > results_view
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename Kokkos::DualView< const scalar_type *, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged >::t_dev diag_view_type
typename matrix_type::local_matrix_device_type local_matrix_type
Kokkos::View< const bool *, memory_space > boundary_nodes_view
typename local_matrix_type::ordinal_type local_ordinal_type
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
DropFunctor(matrix_type &A_, magnitudeType threshold, DistanceFunctorType &dist2_, results_view &results_)
typename local_matrix_type::value_type scalar_type
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
typename matrix_type::local_matrix_device_type local_matrix_type
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
typename local_matrix_type::ordinal_type local_ordinal_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
KokkosKernels::ArithTraits< magnitudeType > mATS
typename local_matrix_type::memory_space memory_space
Kokkos::View< local_ordinal_type *, memory_space > block_indices_view_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
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Kokkos::View< const bool *, memory_space > boundary_nodes_view
typename matrix_type::local_matrix_device_type local_matrix_type
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
VectorDropFunctor(matrix_type &A_, matrix_type &mergedA_, magnitudeType threshold, DistanceFunctorType &dist2_, results_view &results_, block_indices_view_type point_to_block_, block_indices_view_type ghosted_point_to_block_)
KokkosKernels::ArithTraits< impl_scalar_type > implATS
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
KOKKOS_FORCEINLINE_FUNCTION magnitudeType distance2(const local_ordinal_type row, const local_ordinal_type col) const
WeightedDistanceFunctor(matrix_type &A, Teuchos::RCP< coords_type > &coords_, weight_type weight_)
typename matrix_type::local_matrix_device_type local_matrix_type
Xpetra::MultiVector< magnitudeType, LocalOrdinal, GlobalOrdinal, Node > coords_type
typename coords_type::dual_view_type_const::t_dev local_coords_type
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getDiagonal(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, DistanceFunctorType &distFunctor)
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getMaxMinusOffDiagonal(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, DistanceFunctorType &distFunctor)
auto make_vector_drop_functor(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A_, Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &mergedA_, typename VectorDropFunctor< Scalar, LocalOrdinal, GlobalOrdinal, Node, DistanceFunctorType, measure >::magnitudeType threshold, DistanceFunctorType &dist2_, typename VectorDropFunctor< Scalar, LocalOrdinal, GlobalOrdinal, Node, DistanceFunctorType, measure >::results_view &results_, typename VectorDropFunctor< Scalar, LocalOrdinal, GlobalOrdinal, Node, DistanceFunctorType, measure >::block_indices_view_type point_to_block_, typename VectorDropFunctor< Scalar, LocalOrdinal, GlobalOrdinal, Node, DistanceFunctorType, measure >::block_indices_view_type ghosted_point_to_block_)
auto make_drop_functor(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A_, typename DropFunctor< Scalar, LocalOrdinal, GlobalOrdinal, Node, DistanceFunctorType, measure >::magnitudeType threshold, DistanceFunctorType &dist2_, typename DropFunctor< Scalar, LocalOrdinal, GlobalOrdinal, Node, DistanceFunctorType, measure >::results_view &results_)