MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_VectorDroppingDistanceLaplacian_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_VECTORDROPPINGDISTANCELAPLACIAN_DEF_HPP
11 #define MUELU_VECTORDROPPINGDISTANCELAPLACIAN_DEF_HPP
12 
14 
15 namespace MueLu {
16 
17 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, Misc::StrengthMeasure SoC>
19  matrix_type& mergedA,
20  LocalOrdinal blkPartSize,
21  block_indices_view_type& rowTranslation,
22  block_indices_view_type& colTranslation,
23  results_view& results,
24  rowptr_type& filtered_rowptr,
25  rowptr_type& graph_rowptr,
26  nnz_count_type& nnz,
27  boundary_nodes_type& boundaryNodes,
28  const std::string& droppingMethod,
29  const magnitudeType threshold,
30  const bool aggregationMayCreateDirichlet,
31  const bool symmetrizeDroppedGraph,
32  const bool useBlocking,
33  const std::string& distanceLaplacianMetric,
34  Teuchos::Array<double>& dlap_weights,
35  LocalOrdinal interleaved_blocksize,
36  Level& level,
37  const Factory& factory) {
39  auto coords = level.template Get<Teuchos::RCP<doubleMultiVector>>("Coordinates", factory.GetFactory("Coordinates").get());
40  if (distanceLaplacianMetric == "unweighted") {
41  auto dist2 = DistanceLaplacian::UnweightedDistanceFunctor(mergedA, coords);
42  runDroppingFunctors_on_dlap_inner(A, mergedA, blkPartSize, rowTranslation, colTranslation, results, filtered_rowptr, graph_rowptr, nnz, boundaryNodes, droppingMethod, threshold, aggregationMayCreateDirichlet, symmetrizeDroppedGraph, useBlocking, dist2, level, factory);
43  } else if (distanceLaplacianMetric == "weighted") {
44  auto k_dlap_weights_host = Kokkos::View<double*, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>>(&dlap_weights[0], dlap_weights.size());
45  auto k_dlap_weights = Kokkos::View<double*>("dlap_weights", k_dlap_weights_host.extent(0));
46  Kokkos::deep_copy(k_dlap_weights, k_dlap_weights_host);
47  auto dist2 = DistanceLaplacian::WeightedDistanceFunctor(mergedA, coords, k_dlap_weights);
48  runDroppingFunctors_on_dlap_inner(A, mergedA, blkPartSize, rowTranslation, colTranslation, results, filtered_rowptr, graph_rowptr, nnz, boundaryNodes, droppingMethod, threshold, aggregationMayCreateDirichlet, symmetrizeDroppedGraph, useBlocking, dist2, level, factory);
49  } else if (distanceLaplacianMetric == "block weighted") {
50  auto k_dlap_weights_host = Kokkos::View<double*, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>>(&dlap_weights[0], dlap_weights.size());
51  auto k_dlap_weights = Kokkos::View<double*>("dlap_weights", k_dlap_weights_host.extent(0));
52  Kokkos::deep_copy(k_dlap_weights, k_dlap_weights_host);
53  auto dist2 = DistanceLaplacian::BlockWeightedDistanceFunctor(mergedA, coords, k_dlap_weights, interleaved_blocksize);
54  runDroppingFunctors_on_dlap_inner(A, mergedA, blkPartSize, rowTranslation, colTranslation, results, filtered_rowptr, graph_rowptr, nnz, boundaryNodes, droppingMethod, threshold, aggregationMayCreateDirichlet, symmetrizeDroppedGraph, useBlocking, dist2, level, factory);
55  } else if (distanceLaplacianMetric == "material") {
56  auto material = level.template Get<Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>>("Material", factory.GetFactory("Material").get());
57 
58  if (factory.IsPrint(Runtime0)) {
59  auto spatialDim = coords->getNumVectors();
60  if (material->getNumVectors() == 1) {
61  factory.GetOStream(Runtime0) << "material scalar mean = " << material->getVector(0)->meanValue() << std::endl;
62  } else {
63  TEUCHOS_TEST_FOR_EXCEPTION(spatialDim * spatialDim != material->getNumVectors(), Exceptions::RuntimeError, "Need \"Material\" to have spatialDim^2 vectors.");
64  {
65  Teuchos::Array<Scalar> means(material->getNumVectors());
66  material->meanValue(means());
67  std::stringstream ss;
68  ss << "material tensor mean =" << std::endl;
69  size_t k = 0;
70  for (size_t i = 0; i < spatialDim; ++i) {
71  ss << " ";
72  for (size_t j = 0; j < spatialDim; ++j) {
73  ss << means[k] << " ";
74  ++k;
75  }
76  ss << std::endl;
77  }
78  factory.GetOStream(Runtime0) << ss.str();
79  }
80  }
81  }
82 
83  if (material->getNumVectors() == 1) {
84  auto dist2 = DistanceLaplacian::ScalarMaterialDistanceFunctor(mergedA, coords, material);
85  runDroppingFunctors_on_dlap_inner(A, mergedA, blkPartSize, rowTranslation, colTranslation, results, filtered_rowptr, graph_rowptr, nnz, boundaryNodes, droppingMethod, threshold, aggregationMayCreateDirichlet, symmetrizeDroppedGraph, useBlocking, dist2, level, factory);
86  } else {
87  auto dist2 = DistanceLaplacian::TensorMaterialDistanceFunctor(mergedA, coords, material);
88  runDroppingFunctors_on_dlap_inner(A, mergedA, blkPartSize, rowTranslation, colTranslation, results, filtered_rowptr, graph_rowptr, nnz, boundaryNodes, droppingMethod, threshold, aggregationMayCreateDirichlet, symmetrizeDroppedGraph, useBlocking, dist2, level, factory);
89  }
90  }
91 }
92 } // namespace MueLu
93 
94 #define MUELU_ETI_GROUP(SC, LO, GO, NO) \
95  MUELU_ETI_SLGN_SoC(MueLu::VectorDroppingDistanceLaplacian, SC, LO, GO, NO)
96 
97 #endif
MueLu::DefaultLocalOrdinal LocalOrdinal
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
typename Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
One-liner description of what is happening.
Kokkos::View< DecisionType *, memory_space > results_view
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
typename local_graph_type::row_map_type::non_const_type rowptr_type
Kokkos::pair< LocalOrdinal, LocalOrdinal > nnz_count_type
static void runDroppingFunctors_on_dlap(matrix_type &A, matrix_type &mergedA, LocalOrdinal blkPartSize, block_indices_view_type &rowTranslation, block_indices_view_type &colTranslation, results_view &results, rowptr_type &filtered_rowptr, rowptr_type &graph_rowptr, nnz_count_type &nnz, boundary_nodes_type &boundaryNodes, const std::string &droppingMethod, const magnitudeType threshold, const bool aggregationMayCreateDirichlet, const bool symmetrizeDroppedGraph, const bool useBlocking, const std::string &distanceLaplacianMetric, Teuchos::Array< double > &dlap_weights, LocalOrdinal interleaved_blocksize, Level &level, const Factory &factory)
size_type size() const
Exception throws to report errors in the internal logical of the program.
typename Kokkos::View< LocalOrdinal *, typename Node::device_type > block_indices_view_type
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Default implementation of FactoryAcceptor::GetFactory()
typename MueLu::LWGraph_kokkos< LocalOrdinal, GlobalOrdinal, Node >::boundary_nodes_type boundary_nodes_type