Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_SingletonFilter_def.hpp
1/*@HEADER
2// ***********************************************************************
3//
4// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
5// Copyright (2009) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38//
39// ***********************************************************************
40//@HEADER
41*/
42
43#ifndef IFPACK2_SINGLETONFILTER_DEF_HPP
44#define IFPACK2_SINGLETONFILTER_DEF_HPP
45#include "Ifpack2_SingletonFilter_decl.hpp"
46
47#include "Tpetra_ConfigDefs.hpp"
48#include "Tpetra_RowMatrix.hpp"
49#include "Tpetra_Map.hpp"
50#include "Tpetra_MultiVector.hpp"
51#include "Tpetra_Vector.hpp"
52
53namespace Ifpack2 {
54
55template<class MatrixType>
56SingletonFilter<MatrixType>::SingletonFilter(const Teuchos::RCP<const Tpetra::RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >& Matrix):
57 A_(Matrix),
58 NumSingletons_(0),
59 NumRows_(0),
60 NumNonzeros_(0),
61 MaxNumEntries_(0),
62 MaxNumEntriesA_(0)
63{
64
65 // use this filter only on serial matrices
66 if (A_->getComm()->getSize() != 1 || A_->getLocalNumRows() != A_->getGlobalNumRows()) {
67 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.");
68 }
69
70 // Number of rows in A
71 size_t NumRowsA_ = A_->getLocalNumRows();
72
73 // tentative value for MaxNumEntries. This is the number of
74 // nonzeros in the local matrix
75 MaxNumEntriesA_ = A_->getLocalMaxNumRowEntries();
76
77 // ExtractMyRowCopy() will use these vectors
78 Kokkos::resize(Indices_,MaxNumEntriesA_);
79 Kokkos::resize(Values_,MaxNumEntriesA_);
80
81 // Initialize reordering vector to -1
82 Reorder_.resize(NumRowsA_);
83 Reorder_.assign(Reorder_.size(),-1);
84
85 // first check how may singletons I do have
86 NumRows_=0;
87 for (size_t i = 0 ; i < NumRowsA_ ; ++i) {
88 size_t Nnz;
89 A_->getLocalRowCopy(i,Indices_,Values_,Nnz);
90 if (Nnz != 1) {
91 Reorder_[i] = NumRows_++;
92 }
93 else {
94 NumSingletons_++;
95 }
96 }
97
98 // build the inverse reordering
99 InvReorder_.resize(NumRows_);
100 for (size_t i = 0 ; i < NumRowsA_ ; ++i) {
101 if (Reorder_[i] < 0)
102 continue;
103 InvReorder_[Reorder_[i]] = i;
104 }
105 NumEntries_.resize(NumRows_);
106 SingletonIndex_.resize(NumSingletons_);
107
108
109 // now compute the nonzeros per row
110 size_t count = 0;
111 for (size_t i = 0 ; i < NumRowsA_ ; ++i) {
112 size_t Nnz;
113 A_->getLocalRowCopy(i,Indices_,Values_,Nnz);
114 LocalOrdinal ii = Reorder_[i];
115 if (ii >= 0) {
116 NumEntries_[ii] = Nnz;
117 NumNonzeros_ += Nnz;
118 if (Nnz > MaxNumEntries_)
119 MaxNumEntries_ = Nnz;
120 }
121 else {
122 SingletonIndex_[count] = i;
123 count++;
124 }
125 }
126
127 // Build the reduced map. This map should be serial
128 ReducedMap_ = Teuchos::rcp( new Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node>(NumRows_,0,A_->getComm()) );
129
130 // and finish up with the diagonal entry
131 Diagonal_ = Teuchos::rcp( new Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>(ReducedMap_) );
132
133 Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> DiagonalA(A_->getRowMap());
134 A_->getLocalDiagCopy(DiagonalA);
135 const Teuchos::ArrayRCP<const Scalar> & DiagonalAview = DiagonalA.get1dView();
136 for (size_t i = 0 ; i < NumRows_ ; ++i) {
137 LocalOrdinal ii = InvReorder_[i];
138 Diagonal_->replaceLocalValue((LocalOrdinal)i,DiagonalAview[ii]);
139 }
140}
141
142template<class MatrixType>
144
145template<class MatrixType>
146Teuchos::RCP<const Teuchos::Comm<int> >
148{
149 return A_->getComm();
150}
151
152
153template<class MatrixType>
154Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
155 typename MatrixType::global_ordinal_type,
156 typename MatrixType::node_type> >
158{
159 return ReducedMap_;
160}
161
162template<class MatrixType>
163Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
164 typename MatrixType::global_ordinal_type,
165 typename MatrixType::node_type> >
167{
168 return ReducedMap_;
169}
170
171template<class MatrixType>
172Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
173 typename MatrixType::global_ordinal_type,
174 typename MatrixType::node_type> >
176{
177 return ReducedMap_;
178}
179
180template<class MatrixType>
181Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type,
182 typename MatrixType::global_ordinal_type,
183 typename MatrixType::node_type> >
185{
186 return ReducedMap_;
187}
188
189template<class MatrixType>
190Teuchos::RCP<const Tpetra::RowGraph<typename MatrixType::local_ordinal_type,
191 typename MatrixType::global_ordinal_type,
192 typename MatrixType::node_type> >
194{
195 throw std::runtime_error("Ifpack2::SingletonFilter: does not support getGraph.");
196}
197
198template<class MatrixType>
200{
201 return NumRows_;
202}
203
204template<class MatrixType>
206{
207 return NumRows_;
208}
209
210template<class MatrixType>
212{
213 return NumRows_;
214}
215
216template<class MatrixType>
218{
219 return NumRows_;
220}
221
222template<class MatrixType>
223typename MatrixType::global_ordinal_type SingletonFilter<MatrixType>::getIndexBase() const
224{
225 return A_->getIndexBase();
226}
227
228template<class MatrixType>
230{
231 return NumNonzeros_;
232}
233
234template<class MatrixType>
236{
237 return NumNonzeros_;
238}
239
240template<class MatrixType>
241size_t SingletonFilter<MatrixType>::getNumEntriesInGlobalRow(GlobalOrdinal /* globalRow */) const
242{
243 throw std::runtime_error("Ifpack2::SingletonFilter does not implement getNumEntriesInGlobalRow.");
244}
245
246template<class MatrixType>
248{
249 return NumEntries_[localRow];
250}
251
252template<class MatrixType>
254{
255 return MaxNumEntries_;
256}
257
258template<class MatrixType>
260{
261 return MaxNumEntries_;
262}
263
264template<class MatrixType>
266{
267 return true;
268}
269
270template<class MatrixType>
272{
273 return A_->isLocallyIndexed();
274}
275
276template<class MatrixType>
278{
279 return A_->isGloballyIndexed();
280}
281
282template<class MatrixType>
284{
285 return A_->isFillComplete();
286}
287
288template<class MatrixType>
290getGlobalRowCopy (GlobalOrdinal /*LocalRow*/,
291 nonconst_global_inds_host_view_type &/*Indices*/,
292 nonconst_values_host_view_type &/*Values*/,
293 size_t& /*NumEntries*/) const
294{
295 throw std::runtime_error("Ifpack2::SingletonFilter does not implement getGlobalRowCopy.");
296}
297
298
299template<class MatrixType>
301 getLocalRowCopy (LocalOrdinal LocalRow,
302 nonconst_local_inds_host_view_type &Indices,
303 nonconst_values_host_view_type &Values,
304 size_t& NumEntries) const
305{
306 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.");
307
308 size_t Nnz;
309 LocalOrdinal ARow = InvReorder_[LocalRow];
310 A_->getLocalRowCopy(ARow,Indices_,Values_,Nnz);
311
312 // populate the user's vectors
313 NumEntries = 0;
314 for (size_t i = 0 ; i < Nnz ; ++i) {
315 LocalOrdinal ii = Reorder_[Indices_[i]];
316 if ( ii >= 0) {
317 Indices[NumEntries] = ii;
318 Values[NumEntries] = Values_[i];
319 NumEntries++;
320 }
321 }
322}
323
324
325template<class MatrixType>
326void SingletonFilter<MatrixType>::getGlobalRowView(GlobalOrdinal /* GlobalRow */,
327 global_inds_host_view_type &/*indices*/,
328 values_host_view_type &/*values*/) const
329{
330 throw std::runtime_error("Ifpack2::SingletonFilter: does not support getGlobalRowView.");
331}
332
333
334template<class MatrixType>
335void SingletonFilter<MatrixType>::getLocalRowView(LocalOrdinal /* LocalRow */,
336 local_inds_host_view_type & /*indices*/,
337 values_host_view_type & /*values*/) const
338{
339 throw std::runtime_error("Ifpack2::SingletonFilter: does not support getLocalRowView.");
340}
341
342
343template<class MatrixType>
344void SingletonFilter<MatrixType>::getLocalDiagCopy(Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &diag) const
345{
346 // This is somewhat dubious as to how the maps match.
347 return A_->getLocalDiagCopy(diag);
348}
349
350template<class MatrixType>
351void SingletonFilter<MatrixType>::leftScale(const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& /* x */)
352{
353 throw std::runtime_error("Ifpack2::SingletonFilter does not support leftScale.");
354}
355
356template<class MatrixType>
357void SingletonFilter<MatrixType>::rightScale(const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& /* x */)
358{
359 throw std::runtime_error("Ifpack2::SingletonFilter does not support rightScale.");
360}
361
362template<class MatrixType>
363void SingletonFilter<MatrixType>::apply(const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &X,
364 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y,
365 Teuchos::ETransp mode,
366 Scalar /* alpha */,
367 Scalar /* beta */) const
368{
369 typedef Scalar DomainScalar;
370 typedef Scalar RangeScalar;
371
372 // Note: This isn't AztecOO compliant. But neither was Ifpack's version.
373
374 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != Y.getNumVectors(), std::runtime_error,
375 "Ifpack2::SingletonFilter::apply ERROR: X.getNumVectors() != Y.getNumVectors().");
376
377 RangeScalar zero = Teuchos::ScalarTraits<RangeScalar>::zero();
378 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > x_ptr = X.get2dView();
379 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > y_ptr = Y.get2dViewNonConst();
380
381 Y.putScalar(zero);
382 size_t NumVectors = Y.getNumVectors();
383
384
385 for (size_t i = 0 ; i < NumRows_ ; ++i) {
386 size_t Nnz;
387 // Use this class's getrow to make the below code simpler
388 getLocalRowCopy(i,Indices_,Values_,Nnz);
389 if (mode==Teuchos::NO_TRANS){
390 for (size_t j = 0 ; j < Nnz ; ++j)
391 for (size_t k = 0 ; k < NumVectors ; ++k)
392 y_ptr[k][i] += (RangeScalar)Values_[j] * (RangeScalar)x_ptr[k][Indices_[j]];
393 }
394 else if (mode==Teuchos::TRANS){
395 for (size_t j = 0 ; j < Nnz ; ++j)
396 for (size_t k = 0 ; k < NumVectors ; ++k)
397 y_ptr[k][Indices_[j]] += (RangeScalar)Values_[j] * (RangeScalar)x_ptr[k][i];
398 }
399 else { //mode==Teuchos::CONJ_TRANS
400 for (size_t j = 0 ; j < Nnz ; ++j)
401 for (size_t k = 0 ; k < NumVectors ; ++k)
402 y_ptr[k][Indices_[j]] += Teuchos::ScalarTraits<RangeScalar>::conjugate((RangeScalar)Values_[j]) * (RangeScalar)x_ptr[k][i];
403 }
404 }
405}
406
407template<class MatrixType>
409{
410 return true;
411}
412
413template<class MatrixType>
415{
416 return false;
417}
418
419template<class MatrixType>
420void SingletonFilter<MatrixType>::SolveSingletons(const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
421 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
422{
423 this->template SolveSingletonsTempl<Scalar,Scalar>(RHS, LHS);
424}
425
426template<class MatrixType>
427template<class DomainScalar, class RangeScalar>
428void SingletonFilter<MatrixType>::SolveSingletonsTempl(const Tpetra::MultiVector<DomainScalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
429 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
430{
431 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > RHS_ptr = RHS.get2dView();
432 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > LHS_ptr = LHS.get2dViewNonConst();
433
434 for (size_t i = 0 ; i < NumSingletons_ ; ++i) {
435 LocalOrdinal ii = SingletonIndex_[i];
436 // get the diagonal value for the singleton
437 size_t Nnz;
438 A_->getLocalRowCopy(ii,Indices_,Values_,Nnz);
439 for (size_t j = 0 ; j < Nnz ; ++j) {
440 if (Indices_[j] == ii) {
441 for (size_t k = 0 ; k < LHS.getNumVectors() ; ++k)
442 LHS_ptr[k][ii] = (RangeScalar)RHS_ptr[k][ii] / (RangeScalar)Values_[j];
443 }
444 }
445 }
446}
447
448template<class MatrixType>
449void SingletonFilter<MatrixType>::CreateReducedRHS(const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS,
450 const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
451 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS)
452{
453 this->template CreateReducedRHSTempl<Scalar,Scalar>(LHS, RHS, ReducedRHS);
454}
455
456template<class MatrixType>
457template<class DomainScalar, class RangeScalar>
458void SingletonFilter<MatrixType>::CreateReducedRHSTempl(const Tpetra::MultiVector<DomainScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS,
459 const Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
460 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS)
461{
462 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const RangeScalar > > RHS_ptr = RHS.get2dView();
463 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > LHS_ptr = LHS.get2dView();
464 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > ReducedRHS_ptr = ReducedRHS.get2dViewNonConst();
465
466 size_t NumVectors = LHS.getNumVectors();
467
468 for (size_t i = 0 ; i < NumRows_ ; ++i)
469 for (size_t k = 0 ; k < NumVectors ; ++k)
470 ReducedRHS_ptr[k][i] = RHS_ptr[k][InvReorder_[i]];
471
472 for (size_t i = 0 ; i < NumRows_ ; ++i) {
473 LocalOrdinal ii = InvReorder_[i];
474 size_t Nnz;
475 A_->getLocalRowCopy(ii,Indices_,Values_,Nnz);
476
477 for (size_t j = 0 ; j < Nnz ; ++j) {
478 if (Reorder_[Indices_[j]] == -1) {
479 for (size_t k = 0 ; k < NumVectors ; ++k)
480 ReducedRHS_ptr[k][i] -= (RangeScalar)Values_[j] * (RangeScalar)LHS_ptr[k][Indices_[j]];
481 }
482 }
483 }
484}
485
486template<class MatrixType>
487void SingletonFilter<MatrixType>::UpdateLHS(const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedLHS,
488 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
489{
490 this->template UpdateLHSTempl<Scalar,Scalar>(ReducedLHS, LHS);
491}
492
493template<class MatrixType>
494template<class DomainScalar, class RangeScalar>
495void SingletonFilter<MatrixType>::UpdateLHSTempl(const Tpetra::MultiVector<DomainScalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedLHS,
496 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
497{
498
499 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > LHS_ptr = LHS.get2dViewNonConst();
500 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > ReducedLHS_ptr = ReducedLHS.get2dView();
501
502 for (size_t i = 0 ; i < NumRows_ ; ++i)
503 for (size_t k = 0 ; k < LHS.getNumVectors() ; ++k)
504 LHS_ptr[k][InvReorder_[i]] = (RangeScalar)ReducedLHS_ptr[k][i];
505}
506
507template<class MatrixType>
508typename SingletonFilter<MatrixType>::mag_type SingletonFilter<MatrixType>::getFrobeniusNorm() const
509{
510 throw std::runtime_error("Ifpack2::SingletonFilter does not implement getFrobeniusNorm.");
511}
512
513} // namespace Ifpack2
514
515#define IFPACK2_SINGLETONFILTER_INSTANT(S,LO,GO,N) \
516 template class Ifpack2::SingletonFilter< Tpetra::RowMatrix<S, LO, GO, N> >;
517
518#endif
Filter based on matrix entries.
Definition Ifpack2_SingletonFilter_decl.hpp:65
SingletonFilter(const Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Matrix)
Constructor.
Definition Ifpack2_SingletonFilter_def.hpp:56
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:351
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:157
virtual size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:235
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:420
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:344
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:335
virtual size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
Definition Ifpack2_SingletonFilter_def.hpp:259
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:217
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:277
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:247
virtual global_size_t getGlobalNumCols() const
Returns the number of global columns in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:205
virtual mag_type getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
Definition Ifpack2_SingletonFilter_def.hpp:508
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:166
virtual global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:229
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:357
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:290
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:184
virtual bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Definition Ifpack2_SingletonFilter_def.hpp:408
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:301
virtual Teuchos::RCP< const Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const
Returns the RowGraph associated with this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:193
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:241
virtual bool hasColMap() const
Indicates whether this matrix has a well-defined column map.
Definition Ifpack2_SingletonFilter_def.hpp:265
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:487
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:449
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:271
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:175
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:363
virtual global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:199
virtual size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
Definition Ifpack2_SingletonFilter_def.hpp:253
virtual GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this matrix.
Definition Ifpack2_SingletonFilter_def.hpp:223
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:326
virtual bool isFillComplete() const
Returns true if fillComplete() has been called.
Definition Ifpack2_SingletonFilter_def.hpp:283
virtual ~SingletonFilter()
Destructor.
Definition Ifpack2_SingletonFilter_def.hpp:143
virtual size_t getLocalNumRows() const
Returns the number of rows owned on the calling node.
Definition Ifpack2_SingletonFilter_def.hpp:211
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns the communicator.
Definition Ifpack2_SingletonFilter_def.hpp:147
virtual bool supportsRowViews() const
Returns true if RowViews are supported.
Definition Ifpack2_SingletonFilter_def.hpp:414
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:74