10#ifndef MUELU_DISTANCELAPLACIANDROPPING_HPP
11#define MUELU_DISTANCELAPLACIANDROPPING_HPP
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"
32template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
35 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
39 using ATS = Kokkos::ArithTraits<scalar_type>;
41 using implATS = Kokkos::ArithTraits<impl_scalar_type>;
43 using magATS = Kokkos::ArithTraits<magnitudeType>;
44 using coords_type = Xpetra::MultiVector<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>;
56 auto importer = A.getCrsGraph()->getImporter();
57 if (!importer.is_null()) {
58 ghostedCoordsMV = Xpetra::MultiVectorFactory<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(),
coordsMV->getNumVectors());
68 KOKKOS_FORCEINLINE_FUNCTION
72 for (
size_t j = 0; j <
coords.extent(1); ++j) {
80template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
83 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
87 using ATS = Kokkos::ArithTraits<scalar_type>;
89 using implATS = Kokkos::ArithTraits<impl_scalar_type>;
91 using magATS = Kokkos::ArithTraits<magnitudeType>;
92 using coords_type = Xpetra::MultiVector<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>;
94 using material_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
113 auto importer = A.getCrsGraph()->getImporter();
114 if (!importer.is_null()) {
115 ghostedCoordsMV = Xpetra::MultiVectorFactory<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(),
coordsMV->getNumVectors());
120 ghostedMaterialMV = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(),
materialMV->getNumVectors());
133 KOKKOS_INLINE_FUNCTION
143 for (
size_t j = 0; j <
coords.extent(1); ++j) {
151 return Kokkos::max(d_row, d_col);
155template <
class local_ordinal_type,
class material_vector_type,
class material_matrix_type>
162 TensorInversion(material_vector_type& material_vector_, material_matrix_type& material_matrix_)
166 KOKKOS_FORCEINLINE_FUNCTION
173 auto matrix_material = Kokkos::subview(
material_matrix, i, Kokkos::ALL(), Kokkos::ALL());
174 KokkosBatched::SerialLU<KokkosBatched::Algo::LU::Unblocked>::invoke(matrix_material);
178template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
181 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
185 using ATS = Kokkos::ArithTraits<scalar_type>;
187 using implATS = Kokkos::ArithTraits<impl_scalar_type>;
189 using magATS = Kokkos::ArithTraits<magnitudeType>;
190 using coords_type = Xpetra::MultiVector<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>;
192 using material_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
214 auto importer = A.getCrsGraph()->getImporter();
215 if (!importer.is_null()) {
216 ghostedCoordsMV = Xpetra::MultiVectorFactory<magnitudeType, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(),
coordsMV->getNumVectors());
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);
231 ghostedMaterial = material_;
234 using execution_space =
typename Node::execution_space;
235 using range_type = Kokkos::RangePolicy<LocalOrdinal, execution_space>;
238 auto lclMaterial = ghostedMaterial->getDeviceLocalView(Xpetra::Access::ReadOnly);
242 Kokkos::parallel_for(
"MueLu:TensorMaterialDistanceFunctor::inversion", range_type(0, lclMaterial.extent(0)), functor);
246 KOKKOS_INLINE_FUNCTION
255 auto matrix_row_material = Kokkos::subview(
material, row, Kokkos::ALL(), Kokkos::ALL());
256 auto dist = Kokkos::subview(
lcl_dist, row, Kokkos::ALL());
258 for (
size_t j = 0; j <
coords.extent(1); ++j) {
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);
265 for (
size_t j = 0; j <
coords.extent(1); ++j) {
273 auto matrix_col_material = Kokkos::subview(
material, col, Kokkos::ALL(), Kokkos::ALL());
274 auto dist = Kokkos::subview(
lcl_dist, row, Kokkos::ALL());
276 for (
size_t j = 0; j <
coords.extent(1); ++j) {
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);
283 for (
size_t j = 0; j <
coords.extent(1); ++j) {
288 return Kokkos::max(implATS::magnitude(d_row), implATS::magnitude(d_col));
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;
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>;
310 auto diag = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A.getRowMap(), 1);
312 auto lclA = A.getLocalMatrixDevice();
313 auto lclDiag = diag->getDeviceLocalView(Xpetra::Access::OverwriteAll);
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;
323 impl_scalar_type d2 = implATS::zero();
324 for (local_ordinal_type colID = 0; colID < length; colID++) {
325 auto col = rowView.colidx(colID);
327 d = distFunctor.distance2(row, col);
328 d2 += implATS::one() / d;
331 lclDiag(row, 0) = d2;
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);
354template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node,
class DistanceFunctorType>
357 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
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;
367 using ATS = Kokkos::ArithTraits<scalar_type>;
381 :
A(A_.getLocalMatrixDevice())
386 auto lclDiag2d =
diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
387 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
390 KOKKOS_FORCEINLINE_FUNCTION
392 auto row =
A.rowConst(rlid);
393 const size_t offset =
A.graph.row_map(rlid);
395 auto clid = row.colidx(k);
399 val =
one /
dist2.distance2(rlid, clid);
403 auto aiiajj = ATS::magnitude(
diag(rlid)) * ATS::magnitude(
diag(clid));
404 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val);
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
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
DistanceFunctorType dist2
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
Teuchos::RCP< diag_vec_type > diagVec
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
Kokkos::ArithTraits< scalar_type > ATS
DropFunctor(matrix_type &A_, magnitudeType threshold, DistanceFunctorType &dist2_, results_view &results_)
typename ATS::magnitudeType magnitudeType
Teuchos::RCP< material_type > ghostedMaterialMV
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > material_type
Teuchos::RCP< material_type > materialMV
Xpetra::MultiVector< magnitudeType, LocalOrdinal, GlobalOrdinal, Node > coords_type
KOKKOS_INLINE_FUNCTION magnitudeType distance2(const local_ordinal_type row, const local_ordinal_type col) const
Teuchos::RCP< coords_type > coordsMV
typename implATS::magnitudeType magnitudeType
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_)
local_material_type ghostedMaterial
LocalOrdinal local_ordinal_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
local_coords_type ghostedCoords
typename local_matrix_type::value_type scalar_type
typename material_type::dual_view_type_const::t_dev local_material_type
Kokkos::ArithTraits< magnitudeType > magATS
typename matrix_type::local_matrix_type local_matrix_type
Kokkos::ArithTraits< scalar_type > ATS
local_material_type material
Teuchos::RCP< coords_type > ghostedCoordsMV
typename ATS::val_type impl_scalar_type
Kokkos::ArithTraits< impl_scalar_type > implATS
TensorInversion(material_vector_type &material_vector_, material_matrix_type &material_matrix_)
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type i) const
material_matrix_type material_matrix
material_vector_type material_vector
Teuchos::RCP< coords_type > coordsMV
Kokkos::View< impl_scalar_type **, memory_space > local_dist_type
Kokkos::ArithTraits< scalar_type > ATS
typename local_matrix_type::value_type scalar_type
typename coords_type::dual_view_type_const::t_dev local_coords_type
typename ATS::val_type impl_scalar_type
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > material_type
local_coords_type ghostedCoords
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
local_material_type material
Kokkos::View< impl_scalar_type ***, memory_space > local_material_type
Kokkos::ArithTraits< impl_scalar_type > implATS
typename local_matrix_type::memory_space memory_space
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_)
typename matrix_type::local_matrix_type local_matrix_type
Kokkos::ArithTraits< magnitudeType > magATS
Teuchos::RCP< coords_type > ghostedCoordsMV
LocalOrdinal local_ordinal_type
typename implATS::magnitudeType magnitudeType
KOKKOS_FORCEINLINE_FUNCTION magnitudeType distance2(const local_ordinal_type row, const local_ordinal_type col) const
local_coords_type ghostedCoords
Xpetra::MultiVector< magnitudeType, LocalOrdinal, GlobalOrdinal, Node > coords_type
Kokkos::ArithTraits< impl_scalar_type > implATS
Kokkos::ArithTraits< scalar_type > ATS
Teuchos::RCP< coords_type > ghostedCoordsMV
typename local_matrix_type::value_type scalar_type
typename coords_type::dual_view_type_const::t_dev local_coords_type
typename matrix_type::local_matrix_type local_matrix_type
Teuchos::RCP< coords_type > coordsMV
Kokkos::ArithTraits< magnitudeType > magATS
LocalOrdinal local_ordinal_type
typename implATS::magnitudeType magnitudeType
UnweightedDistanceFunctor(matrix_type &A, Teuchos::RCP< coords_type > &coords_)
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
typename ATS::val_type impl_scalar_type
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getDiagonal(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, DistanceFunctorType &distFunctor)