MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_TpetraOperatorAsRowMatrix.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#ifndef MUELU_TPETRAOPERATORASROWMATRIX_HPP
10#define MUELU_TPETRAOPERATORASROWMATRIX_HPP
11
12#include <Teuchos_RCP.hpp>
13#include <Tpetra_RowMatrix.hpp>
14#include <MueLu_Exceptions.hpp>
15
16namespace MueLu {
17
18template <class Scalar = Tpetra::Operator<>::scalar_type,
19 class LocalOrdinal = typename Tpetra::Operator<Scalar>::local_ordinal_type,
20 class GlobalOrdinal = typename Tpetra::Operator<Scalar, LocalOrdinal>::global_ordinal_type,
21 class Node = typename Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
22class TpetraOperatorAsRowMatrix : public Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
23 public:
24 using op_type = Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
25 using vec_type = Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
26
28 using row_matrix_type = Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
29
30 using impl_scalar_type = typename row_matrix_type::impl_scalar_type;
31#if KOKKOS_VERSION >= 40799
32 using mag_type = typename KokkosKernels::ArithTraits<impl_scalar_type>::mag_type;
33#else
34 using mag_type = typename Kokkos::ArithTraits<impl_scalar_type>::mag_type;
35#endif
36
38 typename row_matrix_type::local_inds_device_view_type;
40 typename row_matrix_type::local_inds_host_view_type;
42 typename row_matrix_type::nonconst_local_inds_host_view_type;
43
45 typename row_matrix_type::global_inds_device_view_type;
47 typename row_matrix_type::global_inds_host_view_type;
49 typename row_matrix_type::nonconst_global_inds_host_view_type;
50
52 typename row_matrix_type::values_device_view_type;
54 typename row_matrix_type::values_host_view_type;
56 typename row_matrix_type::nonconst_values_host_view_type;
57
59
60
62 TpetraOperatorAsRowMatrix(const RCP<op_type>& op)
63 : op_(op)
64 , diag_(Teuchos::null) {}
65
66 TpetraOperatorAsRowMatrix(const RCP<op_type>& op,
67 const RCP<vec_type>& diag)
68 : op_(op)
69 , diag_(diag) {}
70
72 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getDomainMap() const {
73 return op_->getDomainMap();
74 }
75
77 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getRangeMap() const {
78 return op_->getRangeMap();
79 }
80
82
86 void apply(const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& X,
87 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Y,
88 Teuchos::ETransp mode = Teuchos::NO_TRANS,
89 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
90 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero()) const {
91 op_->apply(X, Y, mode, alpha, beta);
92 }
93
94 // Fake RowMatrix interface
95 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getRowMap() const {
96 return op_->getRangeMap();
97 }
98
99 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getColMap() const {
100 throw MueLu::Exceptions::RuntimeError("Not implemented.");
101 }
102
103 typename row_matrix_type::local_ordinal_type getBlockSize() const {
104 throw MueLu::Exceptions::RuntimeError("Not implemented.");
105 }
106
107 Teuchos::RCP<const Teuchos::Comm<int> > getComm() const {
108 return op_->getDomainMap()->getComm();
109 }
110
111 Teuchos::RCP<const Tpetra::RowGraph<LocalOrdinal, GlobalOrdinal, Node> > getGraph() const {
112 throw MueLu::Exceptions::RuntimeError("Not implemented.");
113 }
114
115 Tpetra::global_size_t getGlobalNumRows() const {
116 return getRowMap()->getGlobalNumElements();
117 }
118
119 Tpetra::global_size_t getGlobalNumCols() const {
120 return getDomainMap()->getGlobalNumElements();
121 }
122
123 size_t getLocalNumRows() const {
124 return getRowMap()->getLocalNumElements();
125 }
126
127 size_t getLocalNumCols() const {
128 throw MueLu::Exceptions::RuntimeError("Not implemented.");
129 }
130
132 throw MueLu::Exceptions::RuntimeError("Not implemented.");
133 }
134
135 Tpetra::global_size_t getGlobalNumEntries() const {
136 return 0;
137 }
138
139 size_t getLocalNumEntries() const {
140 throw MueLu::Exceptions::RuntimeError("Not implemented.");
141 }
142
143 size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const {
144 throw MueLu::Exceptions::RuntimeError("Not implemented.");
145 }
146
147 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const {
148 throw MueLu::Exceptions::RuntimeError("Not implemented.");
149 }
150
152 throw MueLu::Exceptions::RuntimeError("Not implemented.");
153 }
154
156 throw MueLu::Exceptions::RuntimeError("Not implemented.");
157 }
158
159 bool hasColMap() const {
160 return false;
161 }
162
163 bool isLocallyIndexed() const {
164 return true;
165 }
166
167 bool isGloballyIndexed() const {
168 return true;
169 }
170
171 bool isFillComplete() const {
172 return true;
173 }
174
175 bool supportsRowViews() const {
176 return false;
177 }
178
179 void
183 size_t& NumEntries) const {
184 throw MueLu::Exceptions::RuntimeError("Not implemented.");
185 }
186
187 void
191 size_t& NumEntries) const {
192 throw MueLu::Exceptions::RuntimeError("Not implemented.");
193 }
194
195 void
198 values_host_view_type& values) const {
199 throw MueLu::Exceptions::RuntimeError("Not implemented.");
200 }
201
202 void
205 values_host_view_type& values) const {
206 throw MueLu::Exceptions::RuntimeError("Not implemented.");
207 }
208
209 void getLocalDiagCopy(Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& diag) const {
210 if (diag_.is_null())
211 throw MueLu::Exceptions::RuntimeError("No diagonal available.");
212 else
213 diag = *diag_;
214 }
215
216 void leftScale(const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x) {
217 throw MueLu::Exceptions::RuntimeError("Not implemented.");
218 }
219
220 void rightScale(const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x) {
221 throw MueLu::Exceptions::RuntimeError("Not implemented.");
222 }
223
225 return 0.;
226 }
227
228 // void describe(Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel verbLevel) const {
229 // using std::setw;
230 // using std::endl;
231 // const size_t numRows = nearField_->getRowMap()->getGlobalNumElements();
232 // const size_t nnzNearField = nearField_->getGlobalNumEntries();
233 // const double nnzNearPerRow = Teuchos::as<double>(nnzNearField)/numRows;
234 // const size_t nnzKernelApprox = kernelApproximations_->pointA_->getGlobalNumEntries();
235 // const size_t numClusterPairs = kernelApproximations_->blockA_->getGlobalNumEntries();
236 // const size_t nnzBasis = basisMatrix_->getGlobalNumEntries();
237 // size_t nnzTransfer = 0;
238 // for (size_t i = 0; i<transferMatrices_.size(); i++)
239 // nnzTransfer += transferMatrices_[i]->pointA_->getGlobalNumEntries();
240 // const size_t nnzTotal = nnzNearField+nnzKernelApprox+nnzBasis+nnzTransfer;
241 // const double nnzTotalPerRow = Teuchos::as<double>(nnzTotal)/numRows;
242 // std::ostringstream oss;
243 // oss << std::left;
244 // oss << setw(9) << "rows" << setw(12) << "nnz(near)" << setw(14) << "nnz(near)/row" << setw(12) << "nnz(basis)" << setw(15) << "#cluster pairs" << setw(12)<< "nnz(kernel)" << setw(14) << "nnz(transfer)" << setw(12) << "nnz(total)" << setw(14) << "nnz(total)/row" << endl;
245 // oss << setw(9) << numRows << setw(12) << nnzNearField << setw(14) << nnzNearPerRow << setw(12) << nnzBasis << setw(15) << numClusterPairs << setw(12) << nnzKernelApprox << setw(14) << nnzTransfer << setw(12) << nnzTotal << setw(14) << nnzTotalPerRow << endl;
246 // out << oss.str();
247 // }
248
249 private:
250 RCP<op_type> op_;
251 RCP<vec_type> diag_;
252};
253
254} // namespace MueLu
255
256#endif
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Exception throws to report errors in the internal logical of the program.
Teuchos::RCP< const Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const
void getGlobalRowCopy(GlobalOrdinal GlobalRow, nonconst_global_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const
typename row_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type
void getLocalRowView(LocalOrdinal LocalRow, local_inds_host_view_type &indices, values_host_view_type &values) const
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
typename row_matrix_type::local_inds_host_view_type local_inds_host_view_type
typename row_matrix_type::values_host_view_type values_host_view_type
typename row_matrix_type::global_inds_host_view_type global_inds_host_view_type
void apply(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
Returns in Y the result of a Tpetra::Operator applied to a Tpetra::MultiVector X.
typename row_matrix_type::local_inds_device_view_type local_inds_device_view_type
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
typename row_matrix_type::global_inds_device_view_type global_inds_device_view_type
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > vec_type
void rightScale(const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
typename row_matrix_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type
typename row_matrix_type::impl_scalar_type impl_scalar_type
Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > op_type
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
typename Kokkos::ArithTraits< impl_scalar_type >::mag_type mag_type
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
typename row_matrix_type::values_device_view_type values_device_view_type
typename row_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type
void getLocalRowCopy(LocalOrdinal LocalRow, nonconst_local_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const
Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > row_matrix_type
The RowMatrix representing the base class of CrsMatrix.
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
TpetraOperatorAsRowMatrix(const RCP< op_type > &op, const RCP< vec_type > &diag)
void getGlobalRowView(GlobalOrdinal GlobalRow, global_inds_host_view_type &indices, values_host_view_type &values) const
TpetraOperatorAsRowMatrix(const RCP< op_type > &op)
Constructor.
void leftScale(const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
row_matrix_type::local_ordinal_type getBlockSize() const
void getLocalDiagCopy(Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Namespace for MueLu classes and methods.