MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_BoundaryDetection.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_BOUNDARYDETECTION_HPP
11#define MUELU_BOUNDARYDETECTION_HPP
12
13#include <cstddef>
14#include <type_traits>
15#include "Kokkos_Core.hpp"
16#include "KokkosKernels_ArithTraits.hpp"
17#include "MueLu_LWGraph_kokkos.hpp"
18#include "MueLu_Utilities.hpp"
19#include "Teuchos_RCP.hpp"
20#include "Xpetra_ConfigDefs.hpp"
21#include "Xpetra_CrsGraph.hpp"
22#include "Xpetra_MultiVector.hpp"
23
25
33template <class local_matrix_type>
35 private:
36 using scalar_type = typename local_matrix_type::value_type;
37 using local_ordinal_type = typename local_matrix_type::ordinal_type;
38 using memory_space = typename local_matrix_type::memory_space;
39
40 using ATS = KokkosKernels::ArithTraits<scalar_type>;
41 using magnitudeType = typename ATS::magnitudeType;
42 using boundary_nodes_view = Kokkos::View<bool*, memory_space>;
43
44 local_matrix_type A;
48
49 public:
50 PointDirichletFunctor(local_matrix_type& A_, boundary_nodes_view boundaryNodes_, magnitudeType dirichletThreshold_, local_ordinal_type dirichletNonzeroThreshold_)
51 : A(A_)
52 , boundaryNodes(boundaryNodes_)
53 , dirichletThreshold(dirichletThreshold_)
54 , dirichletNonzeroThreshold(dirichletNonzeroThreshold_) {}
55
56 KOKKOS_FORCEINLINE_FUNCTION
57 void operator()(const local_ordinal_type rlid) const {
58 auto row = A.rowConst(rlid);
59 local_ordinal_type nnz = 0;
60 for (local_ordinal_type k = 0; k < row.length; ++k) {
61 local_ordinal_type clid = row.colidx(k);
62 scalar_type val = row.value(k);
63 if ((rlid != static_cast<local_ordinal_type>(clid)) && (ATS::magnitude(val) > dirichletThreshold)) {
64 ++nnz;
65 if (nnz == dirichletNonzeroThreshold) {
66 return;
67 }
68 }
69 }
70 boundaryNodes(rlid) = true;
71 }
72};
73
83template <class local_matrix_type, bool useGreedyDirichlet>
85 private:
86 using scalar_type = typename local_matrix_type::value_type;
87 using local_ordinal_type = typename local_matrix_type::ordinal_type;
88 using memory_space = typename local_matrix_type::memory_space;
89
90 using ATS = KokkosKernels::ArithTraits<scalar_type>;
91 using magnitudeType = typename ATS::magnitudeType;
92 using boundary_nodes_view = Kokkos::View<bool*, memory_space>;
93
94 local_matrix_type A;
99
100 public:
101 VectorDirichletFunctor(local_matrix_type& A_, local_ordinal_type blockSize_, boundary_nodes_view boundaryNodes_, magnitudeType dirichletThreshold_, local_ordinal_type dirichletNonzeroThreshold_)
102 : A(A_)
103 , blockSize(blockSize_)
104 , boundaryNodes(boundaryNodes_)
105 , dirichletThreshold(dirichletThreshold_)
106 , dirichletNonzeroThreshold(dirichletNonzeroThreshold_) {}
107
108 KOKKOS_FORCEINLINE_FUNCTION
109 void operator()(const local_ordinal_type rblid) const {
110 for (local_ordinal_type rlid = rblid * blockSize; rlid < (rblid + 1) * blockSize; ++rlid) {
111 auto row = A.rowConst(rlid);
112 local_ordinal_type nnz = 0;
113 bool rowIsDirichlet = true;
114 for (local_ordinal_type k = 0; k < row.length; ++k) {
115 auto clid = row.colidx(k);
116 auto val = row.value(k);
117 if ((rlid != static_cast<local_ordinal_type>(clid)) && (ATS::magnitude(val) > dirichletThreshold)) {
118 ++nnz;
119 if (nnz == dirichletNonzeroThreshold) {
120 rowIsDirichlet = false;
121 break;
122 }
123 }
124 }
125 if constexpr (useGreedyDirichlet) {
126 if (rowIsDirichlet) {
127 boundaryNodes(rblid) = true;
128 return;
129 }
130 } else {
131 if (!rowIsDirichlet) {
132 return;
133 }
134 }
135 }
136 if constexpr (!useGreedyDirichlet)
137 boundaryNodes(rblid) = true;
138 }
139};
140
148template <class local_matrix_type>
150 private:
151 using scalar_type = typename local_matrix_type::value_type;
152 using local_ordinal_type = typename local_matrix_type::ordinal_type;
153 using memory_space = typename local_matrix_type::memory_space;
154
155 using ATS = KokkosKernels::ArithTraits<scalar_type>;
156 using magnitudeType = typename ATS::magnitudeType;
157 using magATS = KokkosKernels::ArithTraits<magnitudeType>;
158 using boundary_nodes_view = Kokkos::View<bool*, memory_space>;
159
160 local_matrix_type A;
163
164 public:
165 RowSumFunctor(local_matrix_type& A_, boundary_nodes_view boundaryNodes_, magnitudeType rowSumTol_)
166 : A(A_)
167 , boundaryNodes(boundaryNodes_)
168 , rowSumTol(rowSumTol_) {}
169
170 KOKKOS_FORCEINLINE_FUNCTION
171 void operator()(const local_ordinal_type rlid) const {
172 scalar_type rowsum = ATS::zero();
173 scalar_type diagval = ATS::zero();
174 auto row = A.rowConst(rlid);
175 for (local_ordinal_type k = 0; k < row.length; ++k) {
176 auto clid = row.colidx(k);
177 auto val = row.value(k);
178 if (rlid == static_cast<local_ordinal_type>(clid))
179 diagval = val;
180 rowsum += val;
181 }
182 if (ATS::magnitude(rowsum) > ATS::magnitude(diagval) * rowSumTol) {
183 boundaryNodes(rlid) = true;
184 }
185 }
186};
187
192template <class local_matrix_type, class Functor, class... RemainingFunctors>
194 private:
195 using local_ordinal_type = typename local_matrix_type::ordinal_type;
196
197 Functor functor;
198 BoundaryFunctor<local_matrix_type, RemainingFunctors...> remainingFunctors;
199
200 public:
201 BoundaryFunctor(local_matrix_type& A_, Functor& functor_, RemainingFunctors&... remainingFunctors_)
202 : functor(functor_)
203 , remainingFunctors(A_, remainingFunctors_...) {}
204
205 KOKKOS_FUNCTION void operator()(const local_ordinal_type rlid) const {
206 functor(rlid);
207 remainingFunctors(rlid);
208 }
209};
210
211template <class local_matrix_type, class Functor>
212class BoundaryFunctor<local_matrix_type, Functor> {
213 private:
214 using local_ordinal_type = typename local_matrix_type::ordinal_type;
215
216 local_matrix_type A;
217 Functor functor;
218
219 public:
220 BoundaryFunctor(local_matrix_type& A_, Functor& functor_)
221 : A(A_)
222 , functor(functor_) {}
223
224 KOKKOS_FUNCTION void operator()(const local_ordinal_type rlid) const {
225 functor(rlid);
226 }
227};
228
229} // namespace MueLu::BoundaryDetection
230
231#endif
KOKKOS_FUNCTION void operator()(const local_ordinal_type rlid) const
Functor that serially applies sub-functors to rows.
BoundaryFunctor(local_matrix_type &A_, Functor &functor_, RemainingFunctors &... remainingFunctors_)
BoundaryFunctor< local_matrix_type, RemainingFunctors... > remainingFunctors
typename local_matrix_type::ordinal_type local_ordinal_type
KOKKOS_FUNCTION void operator()(const local_ordinal_type rlid) const
typename local_matrix_type::memory_space memory_space
KokkosKernels::ArithTraits< scalar_type > ATS
Kokkos::View< bool *, memory_space > boundary_nodes_view
typename local_matrix_type::ordinal_type local_ordinal_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(const local_ordinal_type rlid) const
PointDirichletFunctor(local_matrix_type &A_, boundary_nodes_view boundaryNodes_, magnitudeType dirichletThreshold_, local_ordinal_type dirichletNonzeroThreshold_)
typename local_matrix_type::value_type scalar_type
Functor for marking nodes as Dirichlet based on rowsum.
RowSumFunctor(local_matrix_type &A_, boundary_nodes_view boundaryNodes_, magnitudeType rowSumTol_)
typename local_matrix_type::value_type scalar_type
KokkosKernels::ArithTraits< magnitudeType > magATS
typename local_matrix_type::memory_space memory_space
Kokkos::View< bool *, memory_space > boundary_nodes_view
KokkosKernels::ArithTraits< scalar_type > ATS
typename local_matrix_type::ordinal_type local_ordinal_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(const local_ordinal_type rlid) const
Functor for marking nodes as Dirichlet in a block operator.
typename local_matrix_type::ordinal_type local_ordinal_type
typename local_matrix_type::memory_space memory_space
Kokkos::View< bool *, memory_space > boundary_nodes_view
typename local_matrix_type::value_type scalar_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(const local_ordinal_type rblid) const
KokkosKernels::ArithTraits< scalar_type > ATS
VectorDirichletFunctor(local_matrix_type &A_, local_ordinal_type blockSize_, boundary_nodes_view boundaryNodes_, magnitudeType dirichletThreshold_, local_ordinal_type dirichletNonzeroThreshold_)