Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_SingletonFilter_def.hpp
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef IFPACK2_SINGLETONFILTER_DEF_HPP
11#define IFPACK2_SINGLETONFILTER_DEF_HPP
12#include "Ifpack2_SingletonFilter_decl.hpp"
13
14#include "Tpetra_ConfigDefs.hpp"
15#include "Tpetra_RowMatrix.hpp"
16#include "Tpetra_Map.hpp"
17#include "Tpetra_MultiVector.hpp"
18#include "Tpetra_Vector.hpp"
19
20namespace Ifpack2 {
21
22template<class MatrixType>
23SingletonFilter<MatrixType>::SingletonFilter(const Teuchos::RCP<const Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& Matrix):
24 A_(Matrix),
25 NumSingletons_(0),
26 NumRows_(0),
27 NumNonzeros_(0),
28 MaxNumEntries_(0),
29 MaxNumEntriesA_(0)
30{
31
32 // use this filter only on serial matrices
33 if (A_->getComm()->getSize() != 1 || A_->getLocalNumRows() != A_->getGlobalNumRows()) {
34 throw std::runtime_error("Ifpack2::SingeltonFilter can be used with Comm().getSize() == 1 only. This class is a tool for Ifpack2_AdditiveSchwarz, and it is not meant to be used otherwise.");
35 }
36
37 // Number of rows in A
38 size_t NumRowsA_ = A_->getLocalNumRows();
39
40 // tentative value for MaxNumEntries. This is the number of
41 // nonzeros in the local matrix
42 MaxNumEntriesA_ = A_->getLocalMaxNumRowEntries();
43
44 // ExtractMyRowCopy() will use these vectors
45 Kokkos::resize(Indices_,MaxNumEntriesA_);
46 Kokkos::resize(Values_,MaxNumEntriesA_);
47
48 // Initialize reordering vector to -1
49 Reorder_.resize(NumRowsA_);
50 Reorder_.assign(Reorder_.size(),-1);
51
52 // first check how may singletons I do have
53 NumRows_=0;
54 for (size_t i = 0 ; i < NumRowsA_ ; ++i) {
55 size_t Nnz;
56 A_->getLocalRowCopy(i,Indices_,Values_,Nnz);
57 if (Nnz != 1) {
58 Reorder_[i] = NumRows_++;
59 }
60 else {
61 NumSingletons_++;
62 }
63 }
64
65 // build the inverse reordering
66 InvReorder_.resize(NumRows_);
67 for (size_t i = 0 ; i < NumRowsA_ ; ++i) {
68 if (Reorder_[i] < 0)
69 continue;
70 InvReorder_[Reorder_[i]] = i;
71 }
72 NumEntries_.resize(NumRows_);
73 SingletonIndex_.resize(NumSingletons_);
74
75
76 // now compute the nonzeros per row
77 size_t count = 0;
78 for (size_t i = 0 ; i < NumRowsA_ ; ++i) {
79 size_t Nnz;
80 A_->getLocalRowCopy(i,Indices_,Values_,Nnz);
81 LocalOrdinal ii = Reorder_[i];
82 if (ii >= 0) {
83 NumEntries_[ii] = Nnz;
84 NumNonzeros_ += Nnz;
85 if (Nnz > MaxNumEntries_)
86 MaxNumEntries_ = Nnz;
87 }
88 else {
89 SingletonIndex_[count] = i;
90 count++;
91 }
92 }
93
94 // Build the reduced map. This map should be serial
95 ReducedMap_ = Teuchos::rcp( new Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node>(NumRows_,0,A_->getComm()) );
96
97 // and finish up with the diagonal entry
98 Diagonal_ = Teuchos::rcp( new Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>(ReducedMap_) );
99
100 Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> DiagonalA(A_->getRowMap());
101 A_->getLocalDiagCopy(DiagonalA);
102 const Teuchos::ArrayRCP<const Scalar> & DiagonalAview = DiagonalA.get1dView();
103 for (size_t i = 0 ; i < NumRows_ ; ++i) {
104 LocalOrdinal ii = InvReorder_[i];
105 Diagonal_->replaceLocalValue((LocalOrdinal)i,DiagonalAview[ii]);
106 }
107}
108
109template<class MatrixType>
111
112template<class MatrixType>
113Teuchos::RCP<const Teuchos::Comm<int> >
115{
116 return A_->getComm();
117}
118
119
120template<class MatrixType>
121Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
122 typename MatrixType::global_ordinal_type,
123 typename MatrixType::node_type> >
125{
126 return ReducedMap_;
127}
128
129template<class MatrixType>
130Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
131 typename MatrixType::global_ordinal_type,
132 typename MatrixType::node_type> >
134{
135 return ReducedMap_;
136}
137
138template<class MatrixType>
139Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
140 typename MatrixType::global_ordinal_type,
141 typename MatrixType::node_type> >
143{
144 return ReducedMap_;
145}
146
147template<class MatrixType>
148Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
149 typename MatrixType::global_ordinal_type,
150 typename MatrixType::node_type> >
152{
153 return ReducedMap_;
154}
155
156template<class MatrixType>
157Teuchos::RCP<const Tpetra::RowGraph<typename MatrixType::local_ordinal_type,
158 typename MatrixType::global_ordinal_type,
159 typename MatrixType::node_type> >
161{
162 throw std::runtime_error("Ifpack2::SingletonFilter: does not support getGraph.");
163}
164
165template<class MatrixType>
167{
168 return NumRows_;
169}
170
171template<class MatrixType>
173{
174 return NumRows_;
175}
176
177template<class MatrixType>
179{
180 return NumRows_;
181}
182
183template<class MatrixType>
185{
186 return NumRows_;
187}
188
189template<class MatrixType>
190typename MatrixType::global_ordinal_type SingletonFilter<MatrixType>::getIndexBase() const
191{
192 return A_->getIndexBase();
193}
194
195template<class MatrixType>
197{
198 return NumNonzeros_;
199}
200
201template<class MatrixType>
203{
204 return NumNonzeros_;
205}
206
207template<class MatrixType>
208size_t SingletonFilter<MatrixType>::getNumEntriesInGlobalRow(GlobalOrdinal /* globalRow */) const
209{
210 throw std::runtime_error("Ifpack2::SingletonFilter does not implement getNumEntriesInGlobalRow.");
211}
212
213template<class MatrixType>
215{
216 return NumEntries_[localRow];
217}
218
219template<class MatrixType>
221{
222 return MaxNumEntries_;
223}
224
225template<class MatrixType>
227{
228 return MaxNumEntries_;
229}
230
231template<class MatrixType>
232typename MatrixType::local_ordinal_type SingletonFilter<MatrixType>::getBlockSize() const
233{
234 return A_->getBlockSize();
235}
236
237template<class MatrixType>
239{
240 return true;
241}
242
243template<class MatrixType>
245{
246 return A_->isLocallyIndexed();
247}
248
249template<class MatrixType>
251{
252 return A_->isGloballyIndexed();
253}
254
255template<class MatrixType>
257{
258 return A_->isFillComplete();
259}
260
261template<class MatrixType>
263getGlobalRowCopy (GlobalOrdinal /*LocalRow*/,
264 nonconst_global_inds_host_view_type &/*Indices*/,
265 nonconst_values_host_view_type &/*Values*/,
266 size_t& /*NumEntries*/) const
267{
268 throw std::runtime_error("Ifpack2::SingletonFilter does not implement getGlobalRowCopy.");
269}
270
271
272template<class MatrixType>
274 getLocalRowCopy (LocalOrdinal LocalRow,
275 nonconst_local_inds_host_view_type &Indices,
276 nonconst_values_host_view_type &Values,
277 size_t& NumEntries) const
278{
279 TEUCHOS_TEST_FOR_EXCEPTION((LocalRow < 0 || (size_t) LocalRow >= NumRows_ || (size_t) Indices.size() < NumEntries_[LocalRow]), std::runtime_error, "Ifpack2::SingletonFilter::getLocalRowCopy invalid row or array size.");
280
281 size_t Nnz;
282 LocalOrdinal ARow = InvReorder_[LocalRow];
283 A_->getLocalRowCopy(ARow,Indices_,Values_,Nnz);
284
285 // populate the user's vectors
286 NumEntries = 0;
287 for (size_t i = 0 ; i < Nnz ; ++i) {
288 LocalOrdinal ii = Reorder_[Indices_[i]];
289 if ( ii >= 0) {
290 Indices[NumEntries] = ii;
291 Values[NumEntries] = Values_[i];
292 NumEntries++;
293 }
294 }
295}
296
297
298template<class MatrixType>
299void SingletonFilter<MatrixType>::getGlobalRowView(GlobalOrdinal /* GlobalRow */,
300 global_inds_host_view_type &/*indices*/,
301 values_host_view_type &/*values*/) const
302{
303 throw std::runtime_error("Ifpack2::SingletonFilter: does not support getGlobalRowView.");
304}
305
306
307template<class MatrixType>
308void SingletonFilter<MatrixType>::getLocalRowView(LocalOrdinal /* LocalRow */,
309 local_inds_host_view_type & /*indices*/,
310 values_host_view_type & /*values*/) const
311{
312 throw std::runtime_error("Ifpack2::SingletonFilter: does not support getLocalRowView.");
313}
314
315
316template<class MatrixType>
317void SingletonFilter<MatrixType>::getLocalDiagCopy(Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &diag) const
318{
319 // This is somewhat dubious as to how the maps match.
320 return A_->getLocalDiagCopy(diag);
321}
322
323template<class MatrixType>
324void SingletonFilter<MatrixType>::leftScale(const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& /* x */)
325{
326 throw std::runtime_error("Ifpack2::SingletonFilter does not support leftScale.");
327}
328
329template<class MatrixType>
330void SingletonFilter<MatrixType>::rightScale(const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& /* x */)
331{
332 throw std::runtime_error("Ifpack2::SingletonFilter does not support rightScale.");
333}
334
335template<class MatrixType>
336void SingletonFilter<MatrixType>::apply(const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &X,
337 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y,
338 Teuchos::ETransp mode,
339 Scalar /* alpha */,
340 Scalar /* beta */) const
341{
342 typedef Scalar DomainScalar;
343 typedef Scalar RangeScalar;
344
345 // Note: This isn't AztecOO compliant. But neither was Ifpack's version.
346
347 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != Y.getNumVectors(), std::runtime_error,
348 "Ifpack2::SingletonFilter::apply ERROR: X.getNumVectors() != Y.getNumVectors().");
349
350 RangeScalar zero = Teuchos::ScalarTraits<RangeScalar>::zero();
351 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > x_ptr = X.get2dView();
352 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > y_ptr = Y.get2dViewNonConst();
353
354 Y.putScalar(zero);
355 size_t NumVectors = Y.getNumVectors();
356
357
358 for (size_t i = 0 ; i < NumRows_ ; ++i) {
359 size_t Nnz;
360 // Use this class's getrow to make the below code simpler
361 getLocalRowCopy(i,Indices_,Values_,Nnz);
362 if (mode==Teuchos::NO_TRANS){
363 for (size_t j = 0 ; j < Nnz ; ++j)
364 for (size_t k = 0 ; k < NumVectors ; ++k)
365 y_ptr[k][i] += (RangeScalar)Values_[j] * (RangeScalar)x_ptr[k][Indices_[j]];
366 }
367 else if (mode==Teuchos::TRANS){
368 for (size_t j = 0 ; j < Nnz ; ++j)
369 for (size_t k = 0 ; k < NumVectors ; ++k)
370 y_ptr[k][Indices_[j]] += (RangeScalar)Values_[j] * (RangeScalar)x_ptr[k][i];
371 }
372 else { //mode==Teuchos::CONJ_TRANS
373 for (size_t j = 0 ; j < Nnz ; ++j)
374 for (size_t k = 0 ; k < NumVectors ; ++k)
375 y_ptr[k][Indices_[j]] += Teuchos::ScalarTraits<RangeScalar>::conjugate((RangeScalar)Values_[j]) * (RangeScalar)x_ptr[k][i];
376 }
377 }
378}
379
380template<class MatrixType>
382{
383 return true;
384}
385
386template<class MatrixType>
388{
389 return false;
390}
391
392template<class MatrixType>
393void SingletonFilter<MatrixType>::SolveSingletons(const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
394 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
395{
396 this->template SolveSingletonsTempl<Scalar,Scalar>(RHS, LHS);
397}
398
399template<class MatrixType>
400template<class DomainScalar, class RangeScalar>
401void SingletonFilter<MatrixType>::SolveSingletonsTempl(const Tpetra::MultiVector<DomainScalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
402 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
403{
404 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > RHS_ptr = RHS.get2dView();
405 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > LHS_ptr = LHS.get2dViewNonConst();
406
407 for (size_t i = 0 ; i < NumSingletons_ ; ++i) {
408 LocalOrdinal ii = SingletonIndex_[i];
409 // get the diagonal value for the singleton
410 size_t Nnz;
411 A_->getLocalRowCopy(ii,Indices_,Values_,Nnz);
412 for (size_t j = 0 ; j < Nnz ; ++j) {
413 if (Indices_[j] == ii) {
414 for (size_t k = 0 ; k < LHS.getNumVectors() ; ++k)
415 LHS_ptr[k][ii] = (RangeScalar)RHS_ptr[k][ii] / (RangeScalar)Values_[j];
416 }
417 }
418 }
419}
420
421template<class MatrixType>
422void SingletonFilter<MatrixType>::CreateReducedRHS(const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS,
423 const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
424 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS)
425{
426 this->template CreateReducedRHSTempl<Scalar,Scalar>(LHS, RHS, ReducedRHS);
427}
428
429template<class MatrixType>
430template<class DomainScalar, class RangeScalar>
431void SingletonFilter<MatrixType>::CreateReducedRHSTempl(const Tpetra::MultiVector<DomainScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS,
432 const Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
433 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS)
434{
435 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const RangeScalar > > RHS_ptr = RHS.get2dView();
436 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > LHS_ptr = LHS.get2dView();
437 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > ReducedRHS_ptr = ReducedRHS.get2dViewNonConst();
438
439 size_t NumVectors = LHS.getNumVectors();
440
441 for (size_t i = 0 ; i < NumRows_ ; ++i)
442 for (size_t k = 0 ; k < NumVectors ; ++k)
443 ReducedRHS_ptr[k][i] = RHS_ptr[k][InvReorder_[i]];
444
445 for (size_t i = 0 ; i < NumRows_ ; ++i) {
446 LocalOrdinal ii = InvReorder_[i];
447 size_t Nnz;
448 A_->getLocalRowCopy(ii,Indices_,Values_,Nnz);
449
450 for (size_t j = 0 ; j < Nnz ; ++j) {
451 if (Reorder_[Indices_[j]] == -1) {
452 for (size_t k = 0 ; k < NumVectors ; ++k)
453 ReducedRHS_ptr[k][i] -= (RangeScalar)Values_[j] * (RangeScalar)LHS_ptr[k][Indices_[j]];
454 }
455 }
456 }
457}
458
459template<class MatrixType>
460void SingletonFilter<MatrixType>::UpdateLHS(const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedLHS,
461 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
462{
463 this->template UpdateLHSTempl<Scalar,Scalar>(ReducedLHS, LHS);
464}
465
466template<class MatrixType>
467template<class DomainScalar, class RangeScalar>
468void SingletonFilter<MatrixType>::UpdateLHSTempl(const Tpetra::MultiVector<DomainScalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedLHS,
469 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
470{
471
472 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > LHS_ptr = LHS.get2dViewNonConst();
473 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > ReducedLHS_ptr = ReducedLHS.get2dView();
474
475 for (size_t i = 0 ; i < NumRows_ ; ++i)
476 for (size_t k = 0 ; k < LHS.getNumVectors() ; ++k)
477 LHS_ptr[k][InvReorder_[i]] = (RangeScalar)ReducedLHS_ptr[k][i];
478}
479
480template<class MatrixType>
481typename SingletonFilter<MatrixType>::mag_type SingletonFilter<MatrixType>::getFrobeniusNorm() const
482{
483 throw std::runtime_error("Ifpack2::SingletonFilter does not implement getFrobeniusNorm.");
484}
485
486} // namespace Ifpack2
487
488#define IFPACK2_SINGLETONFILTER_INSTANT(S,LO,GO,N) \
489 template class Ifpack2::SingletonFilter< Tpetra::RowMatrix<S, LO, GO, N> >;
490
491#endif
SingletonFilter(const Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Matrix)
Constructor.
Definition Ifpack2_SingletonFilter_def.hpp:23
virtual void leftScale(const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Scales the RowMatrix on the left with the Vector x.
Definition Ifpack2_SingletonFilter_def.hpp:324
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:124
virtual size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:202
virtual void SolveSingletons(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &LHS)
Solve the singleton components of the linear system.
Definition Ifpack2_SingletonFilter_def.hpp:393
virtual void getLocalDiagCopy(Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
Definition Ifpack2_SingletonFilter_def.hpp:317
virtual void getLocalRowView(LocalOrdinal LocalRow, local_inds_host_view_type &indices, values_host_view_type &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
Definition Ifpack2_SingletonFilter_def.hpp:308
virtual size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
Definition Ifpack2_SingletonFilter_def.hpp:226
virtual size_t getLocalNumCols() const
Returns the number of columns needed to apply the forward operator on this node, i....
Definition Ifpack2_SingletonFilter_def.hpp:184
virtual bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
Definition Ifpack2_SingletonFilter_def.hpp:250
virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
Definition Ifpack2_SingletonFilter_def.hpp:214
virtual global_size_t getGlobalNumCols() const
Returns the number of global columns in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:172
virtual mag_type getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
Definition Ifpack2_SingletonFilter_def.hpp:481
virtual LocalOrdinal getBlockSize() const
The number of degrees of freedom per mesh point.
Definition Ifpack2_SingletonFilter_def.hpp:232
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:133
virtual global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:196
virtual void rightScale(const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Scales the RowMatrix on the right with the Vector x.
Definition Ifpack2_SingletonFilter_def.hpp:330
virtual void getGlobalRowCopy(GlobalOrdinal GlobalRow, nonconst_global_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage.
Definition Ifpack2_SingletonFilter_def.hpp:263
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map that describes the range distribution in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:151
virtual bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Definition Ifpack2_SingletonFilter_def.hpp:381
virtual void getLocalRowCopy(LocalOrdinal LocalRow, nonconst_local_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const
Extract a list of entries in a specified local row of the graph. Put into storage allocated by callin...
Definition Ifpack2_SingletonFilter_def.hpp:274
virtual Teuchos::RCP< const Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const
Returns the RowGraph associated with this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:160
virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
Definition Ifpack2_SingletonFilter_def.hpp:208
virtual bool hasColMap() const
Indicates whether this matrix has a well-defined column map.
Definition Ifpack2_SingletonFilter_def.hpp:238
virtual void UpdateLHS(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &ReducedLHS, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &LHS)
Updates a full LHS from a reduces LHS.
Definition Ifpack2_SingletonFilter_def.hpp:460
virtual void CreateReducedRHS(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &LHS, const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &ReducedRHS)
Creates a RHS for the reduced singleton-free system.
Definition Ifpack2_SingletonFilter_def.hpp:422
virtual bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
Definition Ifpack2_SingletonFilter_def.hpp:244
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map that describes the domain distribution in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:142
virtual 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
Computes the operator-multivector application.
Definition Ifpack2_SingletonFilter_def.hpp:336
virtual global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:166
virtual size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
Definition Ifpack2_SingletonFilter_def.hpp:220
virtual GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:190
virtual void getGlobalRowView(GlobalOrdinal GlobalRow, global_inds_host_view_type &indices, values_host_view_type &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
Definition Ifpack2_SingletonFilter_def.hpp:299
virtual bool isFillComplete() const
Returns true if fillComplete() has been called.
Definition Ifpack2_SingletonFilter_def.hpp:256
virtual ~SingletonFilter()
Destructor.
Definition Ifpack2_SingletonFilter_def.hpp:110
virtual size_t getLocalNumRows() const
Returns the number of rows owned on the calling node.
Definition Ifpack2_SingletonFilter_def.hpp:178
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns the communicator.
Definition Ifpack2_SingletonFilter_def.hpp:114
virtual bool supportsRowViews() const
Returns true if RowViews are supported.
Definition Ifpack2_SingletonFilter_def.hpp:387
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41