MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Maxwell_Utils_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_MAXWELL_UTILS_DEF_HPP
11#define MUELU_MAXWELL_UTILS_DEF_HPP
12
13#include <sstream>
14
15#include "MueLu_ConfigDefs.hpp"
16
17#include "Xpetra_Map.hpp"
19#include "Xpetra_MatrixUtils.hpp"
20
22
23#include "MueLu_Level.hpp"
24#include "MueLu_ThresholdAFilterFactory.hpp"
25#include "MueLu_Utilities.hpp"
26#include "MueLu_RAPFactory.hpp"
27
28namespace MueLu {
29
30template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
32 RCP<Matrix>& D0_Matrix_,
33 magnitudeType rowSumTol,
34 bool useKokkos_,
35 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCrowsKokkos_,
36 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCcolsKokkos_,
37 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCdomainKokkos_,
38 int& BCedges_,
39 int& BCnodes_,
40 Teuchos::ArrayRCP<bool>& BCrows_,
41 Teuchos::ArrayRCP<bool>& BCcols_,
42 Teuchos::ArrayRCP<bool>& BCdomain_,
43 bool& allEdgesBoundary_,
44 bool& allNodesBoundary_) {
45 // clean rows associated with boundary conditions
46 // Find rows with only 1 or 2 nonzero entries, record them in BCrows_.
47 // BCrows_[i] is true, iff i is a boundary row
48 // BCcols_[i] is true, iff i is a boundary column
49 int BCedgesLocal = 0;
50 int BCnodesLocal = 0;
51 if (useKokkos_) {
52 BCrowsKokkos_ = Utilities::DetectDirichletRows_kokkos(*SM_Matrix_, Teuchos::ScalarTraits<magnitudeType>::eps(), /*count_twos_as_dirichlet=*/true);
53
54 if (rowSumTol > 0.)
55 Utilities::ApplyRowSumCriterion(*SM_Matrix_, rowSumTol, BCrowsKokkos_);
56
57 BCcolsKokkos_ = Kokkos::View<bool*, typename Node::device_type>(Kokkos::ViewAllocateWithoutInitializing("dirichletCols"), D0_Matrix_->getColMap()->getLocalNumElements());
58 BCdomainKokkos_ = Kokkos::View<bool*, typename Node::device_type>(Kokkos::ViewAllocateWithoutInitializing("dirichletDomains"), D0_Matrix_->getDomainMap()->getLocalNumElements());
59 Utilities::DetectDirichletColsAndDomains(*D0_Matrix_, BCrowsKokkos_, BCcolsKokkos_, BCdomainKokkos_);
60
61 auto BCrowsKokkos = BCrowsKokkos_;
62 Kokkos::parallel_reduce(
63 BCrowsKokkos_.size(), KOKKOS_LAMBDA(int i, int& sum) {
64 if (BCrowsKokkos(i))
65 ++sum;
66 },
67 BCedgesLocal);
68
69 auto BCdomainKokkos = BCdomainKokkos_;
70 Kokkos::parallel_reduce(
71 BCdomainKokkos_.size(), KOKKOS_LAMBDA(int i, int& sum) {
72 if (BCdomainKokkos(i))
73 ++sum;
74 },
75 BCnodesLocal);
76 } else {
77 BCrows_ = Teuchos::arcp_const_cast<bool>(Utilities::DetectDirichletRows(*SM_Matrix_, Teuchos::ScalarTraits<magnitudeType>::eps(), /*count_twos_as_dirichlet=*/true));
78
79 if (rowSumTol > 0.)
80 Utilities::ApplyRowSumCriterion(*SM_Matrix_, rowSumTol, BCrows_);
81
82 BCcols_.resize(D0_Matrix_->getColMap()->getLocalNumElements());
83 BCdomain_.resize(D0_Matrix_->getDomainMap()->getLocalNumElements());
84 Utilities::DetectDirichletColsAndDomains(*D0_Matrix_, BCrows_, BCcols_, BCdomain_);
85
86 for (auto it = BCrows_.begin(); it != BCrows_.end(); ++it)
87 if (*it)
88 BCedgesLocal += 1;
89 for (auto it = BCdomain_.begin(); it != BCdomain_.end(); ++it)
90 if (*it)
91 BCnodesLocal += 1;
92 }
93
94 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCedgesLocal, BCedges_);
95 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCnodesLocal, BCnodes_);
96
97 allEdgesBoundary_ = Teuchos::as<Xpetra::global_size_t>(BCedges_) >= D0_Matrix_->getRangeMap()->getGlobalNumElements();
98 allNodesBoundary_ = Teuchos::as<Xpetra::global_size_t>(BCnodes_) >= D0_Matrix_->getDomainMap()->getGlobalNumElements();
99}
100
101template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
103 RCP<Matrix>& D0_Matrix_,
104 magnitudeType rowSumTol,
105 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCrowsKokkos_,
106 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCcolsKokkos_,
107 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCdomainKokkos_,
108 int& BCedges_,
109 int& BCnodes_,
110 bool& allEdgesBoundary_,
111 bool& allNodesBoundary_) {
112 // clean rows associated with boundary conditions
113 // Find rows with only 1 or 2 nonzero entries, record them in BCrows_.
114 // BCrows_[i] is true, iff i is a boundary row
115 // BCcols_[i] is true, iff i is a boundary column
116 int BCedgesLocal = 0;
117 int BCnodesLocal = 0;
118 {
119 BCrowsKokkos_ = Utilities::DetectDirichletRows_kokkos(*SM_Matrix_, Teuchos::ScalarTraits<magnitudeType>::eps(), /*count_twos_as_dirichlet=*/true);
120
121 if (rowSumTol > 0.)
122 Utilities::ApplyRowSumCriterion(*SM_Matrix_, rowSumTol, BCrowsKokkos_);
123
124 BCcolsKokkos_ = Kokkos::View<bool*, typename Node::device_type>(Kokkos::ViewAllocateWithoutInitializing("dirichletCols"), D0_Matrix_->getColMap()->getLocalNumElements());
125 BCdomainKokkos_ = Kokkos::View<bool*, typename Node::device_type>(Kokkos::ViewAllocateWithoutInitializing("dirichletDomains"), D0_Matrix_->getDomainMap()->getLocalNumElements());
126 Utilities::DetectDirichletColsAndDomains(*D0_Matrix_, BCrowsKokkos_, BCcolsKokkos_, BCdomainKokkos_);
127
128 auto BCrowsKokkos = BCrowsKokkos_;
129 Kokkos::parallel_reduce(
130 BCrowsKokkos_.size(), KOKKOS_LAMBDA(int i, int& sum) {
131 if (BCrowsKokkos(i))
132 ++sum;
133 },
134 BCedgesLocal);
135
136 auto BCdomainKokkos = BCdomainKokkos_;
137 Kokkos::parallel_reduce(
138 BCdomainKokkos_.size(), KOKKOS_LAMBDA(int i, int& sum) {
139 if (BCdomainKokkos(i))
140 ++sum;
141 },
142 BCnodesLocal);
143 }
144
145 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCedgesLocal, BCedges_);
146 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCnodesLocal, BCnodes_);
147
148 allEdgesBoundary_ = Teuchos::as<Xpetra::global_size_t>(BCedges_) >= D0_Matrix_->getRangeMap()->getGlobalNumElements();
149 allNodesBoundary_ = Teuchos::as<Xpetra::global_size_t>(BCnodes_) >= D0_Matrix_->getDomainMap()->getGlobalNumElements();
150}
151
152template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
153RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Maxwell_Utils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
154 removeExplicitZeros(const RCP<Matrix>& A,
155 const magnitudeType tolerance,
156 const bool keepDiagonal) {
157 Level fineLevel;
158 fineLevel.SetFactoryManager(null);
159 fineLevel.SetLevelID(0);
160 fineLevel.Set("A", A);
161 fineLevel.setlib(A->getDomainMap()->lib());
162 RCP<ThresholdAFilterFactory> ThreshFact = rcp(new ThresholdAFilterFactory("A", tolerance, keepDiagonal));
163 fineLevel.Request("A", ThreshFact.get());
164 ThreshFact->Build(fineLevel);
165 return fineLevel.Get<RCP<Matrix> >("A", ThreshFact.get());
166}
167
168template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
170 RCP<Matrix>& D0_Matrix_,
171 RCP<Matrix>& SM_Matrix_,
172 RCP<Matrix>& M1_Matrix_,
173 RCP<Matrix>& Ms_Matrix_) {
174 bool defaultFilter = false;
175
176 // Remove zero entries from D0 if necessary.
177 // In the construction of the prolongator we use the graph of the
178 // matrix, so zero entries mess it up.
179 if (parameterList_.get<bool>("refmaxwell: filter D0", true) && D0_Matrix_->getLocalMaxNumRowEntries() > 2) {
180 D0_Matrix_ = removeExplicitZeros(D0_Matrix_, 1e-8, false);
181
182 // If D0 has too many zeros, maybe SM and M1 do as well.
183 defaultFilter = true;
184 }
185
186 if (!M1_Matrix_.is_null() && parameterList_.get<bool>("refmaxwell: filter M1", defaultFilter)) {
187 RCP<Vector> diag = VectorFactory::Build(M1_Matrix_->getRowMap());
188 // find a reasonable absolute value threshold
189 M1_Matrix_->getLocalDiagCopy(*diag);
190 magnitudeType threshold = 1.0e-8 * diag->normInf();
191
192 M1_Matrix_ = removeExplicitZeros(M1_Matrix_, threshold, true);
193 }
194
195 if (!Ms_Matrix_.is_null() && parameterList_.get<bool>("refmaxwell: filter Ms", defaultFilter)) {
196 RCP<Vector> diag = VectorFactory::Build(Ms_Matrix_->getRowMap());
197 // find a reasonable absolute value threshold
198 Ms_Matrix_->getLocalDiagCopy(*diag);
199 magnitudeType threshold = 1.0e-8 * diag->normInf();
200
201 Ms_Matrix_ = removeExplicitZeros(Ms_Matrix_, threshold, true);
202 }
203
204 if (!SM_Matrix_.is_null() && parameterList_.get<bool>("refmaxwell: filter SM", defaultFilter)) {
205 RCP<Vector> diag = VectorFactory::Build(SM_Matrix_->getRowMap());
206 // find a reasonable absolute value threshold
207 SM_Matrix_->getLocalDiagCopy(*diag);
208 magnitudeType threshold = 1.0e-8 * diag->normInf();
209
210 SM_Matrix_ = removeExplicitZeros(SM_Matrix_, threshold, true);
211 }
212}
213
214template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
216 thresholdedAbs(const RCP<Matrix>& A,
217 const magnitudeType threshold) {
218 using ATS = KokkosKernels::ArithTraits<Scalar>;
219 using impl_Scalar = typename ATS::val_type;
220 using impl_ATS = KokkosKernels::ArithTraits<impl_Scalar>;
221 using range_type = Kokkos::RangePolicy<LO, typename NO::execution_space>;
222
223 const impl_Scalar impl_SC_ONE = impl_ATS::one();
224 const impl_Scalar impl_SC_ZERO = impl_ATS::zero();
225
226 {
227 auto lclMat = A->getLocalMatrixDevice();
228 Kokkos::parallel_for(
229 "thresholdedAbs",
230 range_type(0, lclMat.nnz()),
231 KOKKOS_LAMBDA(const size_t i) {
232 if (impl_ATS::magnitude(lclMat.values(i)) > threshold)
233 lclMat.values(i) = impl_SC_ONE;
234 else
235 lclMat.values(i) = impl_SC_ZERO;
236 });
237 }
238}
239
240template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
242 setMatvecParams(Matrix& A, RCP<ParameterList> matvecParams) {
243 RCP<const Import> xpImporter = A.getCrsGraph()->getImporter();
244 if (!xpImporter.is_null())
245 xpImporter->setDistributorParameters(matvecParams);
246 RCP<const Export> xpExporter = A.getCrsGraph()->getExporter();
247 if (!xpExporter.is_null())
248 xpExporter->setDistributorParameters(matvecParams);
249}
250
251template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
252RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
254 PtAPWrapper(const RCP<Matrix>& A, const RCP<Matrix>& P, ParameterList& params, const std::string& label) {
255 Level fineLevel, coarseLevel;
256 fineLevel.SetFactoryManager(null);
257 coarseLevel.SetFactoryManager(null);
258 coarseLevel.SetPreviousLevel(rcpFromRef(fineLevel));
259 fineLevel.SetLevelID(0);
260 coarseLevel.SetLevelID(1);
261 fineLevel.Set("A", A);
262 coarseLevel.Set("P", P);
263 coarseLevel.setlib(A->getDomainMap()->lib());
264 fineLevel.setlib(A->getDomainMap()->lib());
265 coarseLevel.setObjectLabel(label);
266 fineLevel.setObjectLabel(label);
267
268 RCP<RAPFactory> rapFact = rcp(new RAPFactory());
269 ParameterList rapList = *(rapFact->GetValidParameterList());
270 rapList.set("transpose: use implicit", true);
271 rapList.set("rap: fix zero diagonals", params.get<bool>("rap: fix zero diagonals", true));
272 rapList.set("rap: fix zero diagonals threshold", params.get<double>("rap: fix zero diagonals threshold", Teuchos::ScalarTraits<double>::eps()));
273 rapList.set("rap: triple product", params.get<bool>("rap: triple product", false));
274 rapFact->SetParameterList(rapList);
275
276 coarseLevel.Request("A", rapFact.get());
277
278 return coarseLevel.Get<RCP<Matrix> >("A", rapFact.get());
279}
280
281} // namespace MueLu
282
283#define MUELU_MAXWELL_UTILS_SHORT
284#endif // ifdef MUELU_MAXWELL_UTILS_DEF_HPP
#define MueLu_sumAll(rcpComm, in, out)
Class that holds all level-specific information.
void setlib(Xpetra::UnderlyingLib lib2)
void SetLevelID(int levelID)
Set level number.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
void Request(const FactoryBase &factory)
Increment the storage counter for all the inputs of a factory.
void SetPreviousLevel(const RCP< Level > &previousLevel)
void SetFactoryManager(const RCP< const FactoryManagerBase > &factoryManager)
Set default factories (used internally by Hierarchy::SetLevel()).
static void detectBoundaryConditionsSM(RCP< Matrix > &SM_Matrix, RCP< Matrix > &D0_Matrix, magnitudeType rowSumTol, bool useKokkos_, Kokkos::View< bool *, typename Node::device_type::memory_space > &BCrowsKokkos, Kokkos::View< bool *, typename Node::device_type::memory_space > &BCcolsKokkos, Kokkos::View< bool *, typename Node::device_type::memory_space > &BCdomainKokkos, int &BCedges, int &BCnodes, Teuchos::ArrayRCP< bool > &BCrows, Teuchos::ArrayRCP< bool > &BCcols, Teuchos::ArrayRCP< bool > &BCdomain, bool &allEdgesBoundary, bool &allNodesBoundary)
Detect Dirichlet boundary conditions.
static void thresholdedAbs(const RCP< Matrix > &A, const magnitudeType thresholded)
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
static RCP< Matrix > removeExplicitZeros(const RCP< Matrix > &A, const magnitudeType tolerance, const bool keepDiagonal=true)
Remove explicit zeros.
static void setMatvecParams(Matrix &A, RCP< ParameterList > matvecParams)
Sets matvec params on a matrix.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > PtAPWrapper(const RCP< Matrix > &A, const RCP< Matrix > &P, Teuchos::ParameterList &params, const std::string &label)
Performs an P^T AP.
Factory for building coarse matrices.
Factory for building a thresholded operator.
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude &tol=Teuchos::ScalarTraits< Magnitude >::zero(), bool count_twos_as_dirichlet=false)
Detect Dirichlet rows.
static void ApplyRowSumCriterion(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude rowSumTol, Teuchos::ArrayRCP< bool > &dirichletRows)
Apply Rowsum Criterion.
static Kokkos::View< bool *, typename NO::device_type::memory_space > DetectDirichletRows_kokkos(const Matrix &A, const Magnitude &tol=Teuchos::ScalarTraits< typename Teuchos::ScalarTraits< SC >::magnitudeType >::zero(), const bool count_twos_as_dirichlet=false)
Detect Dirichlet rows.
static void DetectDirichletColsAndDomains(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayRCP< bool > &dirichletRows, Teuchos::ArrayRCP< bool > dirichletCols, Teuchos::ArrayRCP< bool > dirichletDomain)
Detects Dirichlet columns & domains from a list of Dirichlet rows.
Namespace for MueLu classes and methods.
magnitude_type tolerance