MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Maxwell_Utils_decl.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_DECL_HPP
11#define MUELU_MAXWELL_UTILS_DECL_HPP
12
13#include "MueLu_ConfigDefs.hpp"
14#include "MueLu_BaseClass.hpp"
15
16#include "Xpetra_Map_fwd.hpp"
17#include "Xpetra_Matrix_fwd.hpp"
18#include "Xpetra_VectorFactory_fwd.hpp"
19
20#include "MueLu_Level_fwd.hpp"
24
25namespace MueLu {
26
32template <class Scalar,
33 class LocalOrdinal,
34 class GlobalOrdinal,
35 class Node>
37#undef MUELU_MAXWELL_UTILS_SHORT
39
40 public:
41 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitudeType;
42 // typedef typename Teuchos::ScalarTraits<Scalar>::coordinateType coordinateType;
43 // typedef typename Xpetra::MultiVector<coordinateType,LO,GO,NO> RealValuedMultiVector;
44
46 static void detectBoundaryConditionsSM(RCP<Matrix>& SM_Matrix,
47 RCP<Matrix>& D0_Matrix,
48 magnitudeType rowSumTol,
49 bool useKokkos_,
50 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCrowsKokkos,
51 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCcolsKokkos,
52 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCdomainKokkos,
53 int& BCedges,
54 int& BCnodes,
55 Teuchos::ArrayRCP<bool>& BCrows,
56 Teuchos::ArrayRCP<bool>& BCcols,
57 Teuchos::ArrayRCP<bool>& BCdomain,
58 bool& allEdgesBoundary,
59 bool& allNodesBoundary);
60
62 static void detectBoundaryConditionsSM(RCP<Matrix>& SM_Matrix,
63 RCP<Matrix>& D0_Matrix,
64 magnitudeType rowSumTol,
65 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCrowsKokkos,
66 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCcolsKokkos,
67 Kokkos::View<bool*, typename Node::device_type::memory_space>& BCdomainKokkos,
68 int& BCedges,
69 int& BCnodes,
70 bool& allEdgesBoundary,
71 bool& allNodesBoundary);
72
74 static RCP<Matrix> removeExplicitZeros(const RCP<Matrix>& A,
75 const magnitudeType tolerance,
76 const bool keepDiagonal = true,
77 const size_t expectedNNZperRow = 0);
78
79 static void removeExplicitZeros(Teuchos::ParameterList& parameterList,
80 RCP<Matrix>& D0_Matrix,
81 RCP<Matrix>& SM_Matrix,
82 RCP<Matrix>& M1_Matrix,
83 RCP<Matrix>& Ms_Matrix);
84
85 static void removeExplicitZeros(Teuchos::ParameterList& parameterList,
86 RCP<Matrix>& D0_Matrix,
87 RCP<Matrix>& SM_Matrix) {
88 RCP<Matrix> dummy;
89 removeExplicitZeros(parameterList, D0_Matrix, SM_Matrix, dummy, dummy);
90 }
91
92 static void thresholdedAbs(const RCP<Matrix>& A,
93 const magnitudeType thresholded);
94
96 static void setMatvecParams(Matrix& A, RCP<ParameterList> matvecParams);
97
99 static RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
100 PtAPWrapper(const RCP<Matrix>& A, const RCP<Matrix>& P, Teuchos::ParameterList& params, const std::string& label);
101};
102
103} // namespace MueLu
104
105#define MUELU_MAXWELL_UTILS_SHORT
106#endif // MUELU_MAXWELL_UTILS_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Utility functions for Maxwell.
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 void removeExplicitZeros(Teuchos::ParameterList &parameterList, RCP< Matrix > &D0_Matrix, RCP< Matrix > &SM_Matrix)
static RCP< Matrix > removeExplicitZeros(const RCP< Matrix > &A, const magnitudeType tolerance, const bool keepDiagonal=true, const size_t expectedNNZperRow=0)
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.
Verbose class for MueLu classes.
Namespace for MueLu classes and methods.