46#ifndef MUELU_MAXWELL_UTILS_DEF_HPP
47#define MUELU_MAXWELL_UTILS_DEF_HPP
53#include "Xpetra_Map.hpp"
54#include "Xpetra_CrsMatrixUtils.hpp"
55#include "Xpetra_MatrixUtils.hpp"
60#include "MueLu_ThresholdAFilterFactory.hpp"
61#include "MueLu_Utilities.hpp"
62#include "MueLu_RAPFactory.hpp"
68 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
70 RCP<Matrix> & D0_Matrix_,
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_,
78 Teuchos::ArrayRCP<bool> & BCrows_,
79 Teuchos::ArrayRCP<bool> & BCcols_,
80 Teuchos::ArrayRCP<bool> & BCdomain_,
81 bool & allEdgesBoundary_,
82 bool & allNodesBoundary_) {
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());
99 auto BCrowsKokkos=BCrowsKokkos_;
100 Kokkos::parallel_reduce(BCrowsKokkos_.size(), KOKKOS_LAMBDA (
int i,
int & sum) {
105 auto BCdomainKokkos = BCdomainKokkos_;
106 Kokkos::parallel_reduce(BCdomainKokkos_.size(), KOKKOS_LAMBDA (
int i,
int & sum) {
107 if (BCdomainKokkos(i))
117 BCcols_.resize(D0_Matrix_->getColMap()->getLocalNumElements());
118 BCdomain_.resize(D0_Matrix_->getDomainMap()->getLocalNumElements());
121 for (
auto it = BCrows_.begin(); it != BCrows_.end(); ++it)
124 for (
auto it = BCdomain_.begin(); it != BCdomain_.end(); ++it)
129 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCedgesLocal, BCedges_);
130 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCnodesLocal, BCnodes_);
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();
138 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
140 RCP<Matrix> & D0_Matrix_,
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_,
147 bool & allEdgesBoundary_,
148 bool & allNodesBoundary_) {
153 int BCedgesLocal = 0;
154 int BCnodesLocal = 0;
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());
164 auto BCrowsKokkos=BCrowsKokkos_;
165 Kokkos::parallel_reduce(BCrowsKokkos_.size(), KOKKOS_LAMBDA (
int i,
int & sum) {
170 auto BCdomainKokkos = BCdomainKokkos_;
171 Kokkos::parallel_reduce(BCdomainKokkos_.size(), KOKKOS_LAMBDA (
int i,
int & sum) {
172 if (BCdomainKokkos(i))
176 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCedgesLocal, BCedges_);
177 MueLu_sumAll(SM_Matrix_->getRowMap()->getComm(), BCnodesLocal, BCnodes_);
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();
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_) {
192 bool defaultFilter =
false;
197 if (parameterList_.get<
bool>(
"refmaxwell: filter D0",
true) && D0_Matrix_->getLocalMaxNumRowEntries()>2) {
201 fineLevel.
Set(
"A",D0_Matrix_);
202 fineLevel.
setlib(D0_Matrix_->getDomainMap()->lib());
205 fineLevel.
Request(
"A",ThreshFact.get());
206 ThreshFact->Build(fineLevel);
207 D0_Matrix_ = fineLevel.
Get< RCP<Matrix> >(
"A",ThreshFact.get());
210 defaultFilter =
true;
213 if (!M1_Matrix_.is_null() && parameterList_.get<
bool>(
"refmaxwell: filter M1", defaultFilter)) {
214 RCP<Vector> diag = VectorFactory::Build(M1_Matrix_->getRowMap());
216 M1_Matrix_->getLocalDiagCopy(*diag);
222 fineLevel.
Set(
"A",M1_Matrix_);
223 fineLevel.
setlib(M1_Matrix_->getDomainMap()->lib());
225 fineLevel.
Request(
"A",ThreshFact.get());
226 ThreshFact->Build(fineLevel);
227 M1_Matrix_ = fineLevel.
Get< RCP<Matrix> >(
"A",ThreshFact.get());
230 if (!Ms_Matrix_.is_null() && parameterList_.get<
bool>(
"refmaxwell: filter Ms", defaultFilter)) {
231 RCP<Vector> diag = VectorFactory::Build(Ms_Matrix_->getRowMap());
233 Ms_Matrix_->getLocalDiagCopy(*diag);
239 fineLevel.
Set(
"A",Ms_Matrix_);
240 fineLevel.
setlib(Ms_Matrix_->getDomainMap()->lib());
242 fineLevel.
Request(
"A",ThreshFact.get());
243 ThreshFact->Build(fineLevel);
244 Ms_Matrix_ = fineLevel.
Get< RCP<Matrix> >(
"A",ThreshFact.get());
247 if (!SM_Matrix_.is_null() && parameterList_.get<
bool>(
"refmaxwell: filter SM", defaultFilter)) {
248 RCP<Vector> diag = VectorFactory::Build(SM_Matrix_->getRowMap());
250 SM_Matrix_->getLocalDiagCopy(*diag);
256 fineLevel.
Set(
"A",SM_Matrix_);
257 fineLevel.
setlib(SM_Matrix_->getDomainMap()->lib());
259 fineLevel.
Request(
"A",ThreshFact.get());
260 ThreshFact->Build(fineLevel);
261 SM_Matrix_ = fineLevel.
Get< RCP<Matrix> >(
"A",ThreshFact.get());
268 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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);
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 ¶ms, std::string & label) {
284 Level fineLevel, coarseLevel;
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);
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);
305 coarseLevel.
Request(
"A", rapFact.get());
307 return coarseLevel.
Get< RCP<Matrix> >(
"A", rapFact.get());
317#define MUELU_MAXWELL_UTILS_SHORT
#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 ¶meterList, 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 ¶ms, 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.