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// ***********************************************************************
4//
5// MueLu: A package for multigrid based preconditioning
6// Copyright 2012 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact
39// Jonathan Hu (jhu@sandia.gov)
40// Andrey Prokopenko (aprokop@sandia.gov)
41// Ray Tuminaro (rstumin@sandia.gov)
42//
43// ***********************************************************************
44//
45// @HEADER
46#ifndef MUELU_MAXWELL_UTILS_DEF_HPP
47#define MUELU_MAXWELL_UTILS_DEF_HPP
48
49#include <sstream>
50
51#include "MueLu_ConfigDefs.hpp"
52
53#include "Xpetra_Map.hpp"
54#include "Xpetra_CrsMatrixUtils.hpp"
55#include "Xpetra_MatrixUtils.hpp"
56
58
59#include "MueLu_Level.hpp"
60#include "MueLu_ThresholdAFilterFactory.hpp"
61#include "MueLu_Utilities.hpp"
62#include "MueLu_RAPFactory.hpp"
63
64
65namespace MueLu {
66
67
68 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
70 RCP<Matrix> & D0_Matrix_,
71 magnitudeType rowSumTol,
72 bool useKokkos_,
73 Kokkos::View<bool*, typename Node::device_type> & BCrowsKokkos_,
74 Kokkos::View<bool*, typename Node::device_type> & BCcolsKokkos_,
75 Kokkos::View<bool*, typename Node::device_type> & BCdomainKokkos_,
76 int & BCedges_,
77 int & BCnodes_,
78 Teuchos::ArrayRCP<bool> & BCrows_,
79 Teuchos::ArrayRCP<bool> & BCcols_,
80 Teuchos::ArrayRCP<bool> & BCdomain_,
81 bool & allEdgesBoundary_,
82 bool & allNodesBoundary_) {
83 // clean rows associated with boundary conditions
84 // Find rows with only 1 or 2 nonzero entries, record them in BCrows_.
85 // BCrows_[i] is true, iff i is a boundary row
86 // BCcols_[i] is true, iff i is a boundary column
87 int BCedgesLocal = 0;
88 int BCnodesLocal = 0;
89 if (useKokkos_) {
90 BCrowsKokkos_ = Utilities::DetectDirichletRows_kokkos(*SM_Matrix_,Teuchos::ScalarTraits<magnitudeType>::eps(),/*count_twos_as_dirichlet=*/true);
91
92 if (rowSumTol > 0.)
93 Utilities::ApplyRowSumCriterion(*SM_Matrix_, rowSumTol, BCrowsKokkos_);
94
95 BCcolsKokkos_ = Kokkos::View<bool*,typename Node::device_type>(Kokkos::ViewAllocateWithoutInitializing("dirichletCols"), D0_Matrix_->getColMap()->getLocalNumElements());
96 BCdomainKokkos_ = Kokkos::View<bool*,typename Node::device_type>(Kokkos::ViewAllocateWithoutInitializing("dirichletDomains"), D0_Matrix_->getDomainMap()->getLocalNumElements());
97 Utilities::DetectDirichletColsAndDomains(*D0_Matrix_,BCrowsKokkos_,BCcolsKokkos_,BCdomainKokkos_);
98
99 auto BCrowsKokkos=BCrowsKokkos_;
100 Kokkos::parallel_reduce(BCrowsKokkos_.size(), KOKKOS_LAMBDA (int i, int & sum) {
101 if (BCrowsKokkos(i))
102 ++sum;
103 }, BCedgesLocal );
104
105 auto BCdomainKokkos = BCdomainKokkos_;
106 Kokkos::parallel_reduce(BCdomainKokkos_.size(), KOKKOS_LAMBDA (int i, int & sum) {
107 if (BCdomainKokkos(i))
108 ++sum;
109 }, BCnodesLocal);
110 } else
111 {
112 BCrows_ = Teuchos::arcp_const_cast<bool>(Utilities::DetectDirichletRows(*SM_Matrix_,Teuchos::ScalarTraits<magnitudeType>::eps(),/*count_twos_as_dirichlet=*/true));
113
114 if (rowSumTol > 0.)
115 Utilities::ApplyRowSumCriterion(*SM_Matrix_, rowSumTol, BCrows_);
116
117 BCcols_.resize(D0_Matrix_->getColMap()->getLocalNumElements());
118 BCdomain_.resize(D0_Matrix_->getDomainMap()->getLocalNumElements());
119 Utilities::DetectDirichletColsAndDomains(*D0_Matrix_,BCrows_,BCcols_,BCdomain_);
120
121 for (auto it = BCrows_.begin(); it != BCrows_.end(); ++it)
122 if (*it)
123 BCedgesLocal += 1;
124 for (auto it = BCdomain_.begin(); it != BCdomain_.end(); ++it)
125 if (*it)
126 BCnodesLocal += 1;
127 }
128
129 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCedgesLocal, BCedges_);
130 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCnodesLocal, BCnodes_);
131
132
133 allEdgesBoundary_ = Teuchos::as<Xpetra::global_size_t>(BCedges_) >= D0_Matrix_->getRangeMap()->getGlobalNumElements();
134 allNodesBoundary_ = Teuchos::as<Xpetra::global_size_t>(BCnodes_) >= D0_Matrix_->getDomainMap()->getGlobalNumElements();
135 }
136
137
138 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
140 RCP<Matrix> & D0_Matrix_,
141 magnitudeType rowSumTol,
142 Kokkos::View<bool*, typename Node::device_type> & BCrowsKokkos_,
143 Kokkos::View<bool*, typename Node::device_type> & BCcolsKokkos_,
144 Kokkos::View<bool*, typename Node::device_type> & BCdomainKokkos_,
145 int & BCedges_,
146 int & BCnodes_,
147 bool & allEdgesBoundary_,
148 bool & allNodesBoundary_) {
149 // clean rows associated with boundary conditions
150 // Find rows with only 1 or 2 nonzero entries, record them in BCrows_.
151 // BCrows_[i] is true, iff i is a boundary row
152 // BCcols_[i] is true, iff i is a boundary column
153 int BCedgesLocal = 0;
154 int BCnodesLocal = 0;
155 BCrowsKokkos_ = Utilities::DetectDirichletRows_kokkos(*SM_Matrix_,Teuchos::ScalarTraits<magnitudeType>::eps(),/*count_twos_as_dirichlet=*/true);
156
157 if (rowSumTol > 0.)
158 Utilities::ApplyRowSumCriterion(*SM_Matrix_, rowSumTol, BCrowsKokkos_);
159
160 BCcolsKokkos_ = Kokkos::View<bool*,typename Node::device_type>(Kokkos::ViewAllocateWithoutInitializing("dirichletCols"), D0_Matrix_->getColMap()->getLocalNumElements());
161 BCdomainKokkos_ = Kokkos::View<bool*,typename Node::device_type>(Kokkos::ViewAllocateWithoutInitializing("dirichletDomains"), D0_Matrix_->getDomainMap()->getLocalNumElements());
162 Utilities::DetectDirichletColsAndDomains(*D0_Matrix_,BCrowsKokkos_,BCcolsKokkos_,BCdomainKokkos_);
163
164 auto BCrowsKokkos=BCrowsKokkos_;
165 Kokkos::parallel_reduce(BCrowsKokkos_.size(), KOKKOS_LAMBDA (int i, int & sum) {
166 if (BCrowsKokkos(i))
167 ++sum;
168 }, BCedgesLocal );
169
170 auto BCdomainKokkos = BCdomainKokkos_;
171 Kokkos::parallel_reduce(BCdomainKokkos_.size(), KOKKOS_LAMBDA (int i, int & sum) {
172 if (BCdomainKokkos(i))
173 ++sum;
174 }, BCnodesLocal);
175
176 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCedgesLocal, BCedges_);
177 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCnodesLocal, BCnodes_);
178
179
180 allEdgesBoundary_ = Teuchos::as<Xpetra::global_size_t>(BCedges_) >= D0_Matrix_->getRangeMap()->getGlobalNumElements();
181 allNodesBoundary_ = Teuchos::as<Xpetra::global_size_t>(BCnodes_) >= D0_Matrix_->getDomainMap()->getGlobalNumElements();
182 }
183
184
185 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
187 RCP<Matrix> & D0_Matrix_,
188 RCP<Matrix> & SM_Matrix_,
189 RCP<Matrix> & M1_Matrix_,
190 RCP<Matrix> & Ms_Matrix_) {
191
192 bool defaultFilter = false;
193
194 // Remove zero entries from D0 if necessary.
195 // In the construction of the prolongator we use the graph of the
196 // matrix, so zero entries mess it up.
197 if (parameterList_.get<bool>("refmaxwell: filter D0", true) && D0_Matrix_->getLocalMaxNumRowEntries()>2) {
198 Level fineLevel;
199 fineLevel.SetFactoryManager(null);
200 fineLevel.SetLevelID(0);
201 fineLevel.Set("A",D0_Matrix_);
202 fineLevel.setlib(D0_Matrix_->getDomainMap()->lib());
203 // We expect D0 to have entries +-1, so any threshold value will do.
204 RCP<ThresholdAFilterFactory> ThreshFact = rcp(new ThresholdAFilterFactory("A",1.0e-8,/*keepDiagonal=*/false,/*expectedNNZperRow=*/2));
205 fineLevel.Request("A",ThreshFact.get());
206 ThreshFact->Build(fineLevel);
207 D0_Matrix_ = fineLevel.Get< RCP<Matrix> >("A",ThreshFact.get());
208
209 // If D0 has too many zeros, maybe SM and M1 do as well.
210 defaultFilter = true;
211 }
212
213 if (!M1_Matrix_.is_null() && parameterList_.get<bool>("refmaxwell: filter M1", defaultFilter)) {
214 RCP<Vector> diag = VectorFactory::Build(M1_Matrix_->getRowMap());
215 // find a reasonable absolute value threshold
216 M1_Matrix_->getLocalDiagCopy(*diag);
217 magnitudeType threshold = 1.0e-8 * diag->normInf();
218
219 Level fineLevel;
220 fineLevel.SetFactoryManager(null);
221 fineLevel.SetLevelID(0);
222 fineLevel.Set("A",M1_Matrix_);
223 fineLevel.setlib(M1_Matrix_->getDomainMap()->lib());
224 RCP<ThresholdAFilterFactory> ThreshFact = rcp(new ThresholdAFilterFactory("A",threshold,/*keepDiagonal=*/true));
225 fineLevel.Request("A",ThreshFact.get());
226 ThreshFact->Build(fineLevel);
227 M1_Matrix_ = fineLevel.Get< RCP<Matrix> >("A",ThreshFact.get());
228 }
229
230 if (!Ms_Matrix_.is_null() && parameterList_.get<bool>("refmaxwell: filter Ms", defaultFilter)) {
231 RCP<Vector> diag = VectorFactory::Build(Ms_Matrix_->getRowMap());
232 // find a reasonable absolute value threshold
233 Ms_Matrix_->getLocalDiagCopy(*diag);
234 magnitudeType threshold = 1.0e-8 * diag->normInf();
235
236 Level fineLevel;
237 fineLevel.SetFactoryManager(null);
238 fineLevel.SetLevelID(0);
239 fineLevel.Set("A",Ms_Matrix_);
240 fineLevel.setlib(Ms_Matrix_->getDomainMap()->lib());
241 RCP<ThresholdAFilterFactory> ThreshFact = rcp(new ThresholdAFilterFactory("A",threshold,/*keepDiagonal=*/true));
242 fineLevel.Request("A",ThreshFact.get());
243 ThreshFact->Build(fineLevel);
244 Ms_Matrix_ = fineLevel.Get< RCP<Matrix> >("A",ThreshFact.get());
245 }
246
247 if (!SM_Matrix_.is_null() && parameterList_.get<bool>("refmaxwell: filter SM", defaultFilter)) {
248 RCP<Vector> diag = VectorFactory::Build(SM_Matrix_->getRowMap());
249 // find a reasonable absolute value threshold
250 SM_Matrix_->getLocalDiagCopy(*diag);
251 magnitudeType threshold = 1.0e-8 * diag->normInf();
252
253 Level fineLevel;
254 fineLevel.SetFactoryManager(null);
255 fineLevel.SetLevelID(0);
256 fineLevel.Set("A",SM_Matrix_);
257 fineLevel.setlib(SM_Matrix_->getDomainMap()->lib());
258 RCP<ThresholdAFilterFactory> ThreshFact = rcp(new ThresholdAFilterFactory("A",threshold,/*keepDiagonal=*/true));
259 fineLevel.Request("A",ThreshFact.get());
260 ThreshFact->Build(fineLevel);
261 SM_Matrix_ = fineLevel.Get< RCP<Matrix> >("A",ThreshFact.get());
262 }
263
264 }
265
266
267
268 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
270 setMatvecParams(Matrix& A, RCP<ParameterList> matvecParams) {
271 RCP<const Import> xpImporter = A.getCrsGraph()->getImporter();
272 if (!xpImporter.is_null())
273 xpImporter->setDistributorParameters(matvecParams);
274 RCP<const Export> xpExporter = A.getCrsGraph()->getExporter();
275 if (!xpExporter.is_null())
276 xpExporter->setDistributorParameters(matvecParams);
277 }
278
279
280 template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
281 RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
283 PtAPWrapper(const RCP<Matrix>& A,const RCP<Matrix>& P, ParameterList &params, std::string & label) {
284 Level fineLevel, coarseLevel;
285 fineLevel.SetFactoryManager(null);
286 coarseLevel.SetFactoryManager(null);
287 coarseLevel.SetPreviousLevel(rcpFromRef(fineLevel));
288 fineLevel.SetLevelID(0);
289 coarseLevel.SetLevelID(1);
290 fineLevel.Set("A",A);
291 coarseLevel.Set("P",P);
292 coarseLevel.setlib(A->getDomainMap()->lib());
293 fineLevel.setlib(A->getDomainMap()->lib());
294 coarseLevel.setObjectLabel(label);
295 fineLevel.setObjectLabel(label);
296
297 RCP<RAPFactory> rapFact = rcp(new RAPFactory());
298 ParameterList rapList = *(rapFact->GetValidParameterList());
299 rapList.set("transpose: use implicit", true);
300 rapList.set("rap: fix zero diagonals", params.get<bool>("rap: fix zero diagonals", true));
301 rapList.set("rap: fix zero diagonals threshold", params.get<double>("rap: fix zero diagonals threshold", Teuchos::ScalarTraits<double>::eps()));
302 rapList.set("rap: triple product", params.get<bool>("rap: triple product", false));
303 rapFact->SetParameterList(rapList);
304
305 coarseLevel.Request("A", rapFact.get());
306
307 return coarseLevel.Get< RCP<Matrix> >("A", rapFact.get());
308 }
309
310
311
312
313
314
315} // namespace
316
317#define MUELU_MAXWELL_UTILS_SHORT
318#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 removeExplicitZeros(Teuchos::ParameterList &parameterList, RCP< Matrix > &D0_Matrix, RCP< Matrix > &SM_Matrix, RCP< Matrix > &M1_Matrix, RCP< Matrix > &Ms_Matrix)
Remove explicit zeros.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > PtAPWrapper(const RCP< Matrix > &A, const RCP< Matrix > &P, Teuchos::ParameterList &params, std::string &label)
Performs an P^T AP.
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitudeType
static void detectBoundaryConditionsSM(RCP< Matrix > &SM_Matrix, RCP< Matrix > &D0_Matrix, magnitudeType rowSumTol, bool useKokkos_, Kokkos::View< bool *, typename Node::device_type > &BCrowsKokkos, Kokkos::View< bool *, typename Node::device_type > &BCcolsKokkos, Kokkos::View< bool *, typename Node::device_type > &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 setMatvecParams(Matrix &A, RCP< ParameterList > matvecParams)
Sets matvec params on a matrix.
Factory for building coarse matrices.
Factory for building a thresholded operator.
static Kokkos::View< bool *, typename NO::device_type > DetectDirichletRows_kokkos(const Matrix &A, const Magnitude &tol=Teuchos::ScalarTraits< typename Teuchos::ScalarTraits< SC >::magnitudeType >::zero(), const bool count_twos_as_dirichlet=false)
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &A, const Magnitude &tol=Teuchos::ScalarTraits< Magnitude >::zero(), bool count_twos_as_dirichlet=false)
static void ApplyRowSumCriterion(const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &A, const Magnitude rowSumTol, Teuchos::ArrayRCP< bool > &dirichletRows)
static void DetectDirichletColsAndDomains(const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &A, const Teuchos::ArrayRCP< bool > &dirichletRows, Teuchos::ArrayRCP< bool > dirichletCols, Teuchos::ArrayRCP< bool > dirichletDomain)
Namespace for MueLu classes and methods.