Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_ReorderFilter_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_REORDERFILTER_DEF_HPP
44 #define IFPACK2_REORDERFILTER_DEF_HPP
45 #include "Ifpack2_ReorderFilter_decl.hpp"
46 #include <vector>
47 
48 #include "Tpetra_ConfigDefs.hpp"
49 #include "Tpetra_RowMatrix.hpp"
50 #include "Tpetra_Map.hpp"
51 #include "Tpetra_MultiVector.hpp"
52 #include "Tpetra_Vector.hpp"
53 
54 namespace Ifpack2 {
55 
56 template<class MatrixType>
60  const Teuchos::ArrayRCP<local_ordinal_type>& reverseperm)
61  : A_ (A),
62  perm_ (perm),
63  reverseperm_ (reverseperm)
64 {
66  A_.is_null (), std::invalid_argument,
67  "Ifpack2::ReorderFilter: The input matrix is null.");
68 
69  // use this filter only on serial matrices
71  A_->getComm()->getSize() != 1, std::invalid_argument,
72  "Ifpack2::ReorderFilter: This class may only be used if the input matrix's "
73  "communicator has one process. This class is an implementation detail of "
74  "Ifpack2::AdditiveSchwarz, and it is not meant to be used otherwise.");
75 
77  A_->getLocalNumRows () != A_->getGlobalNumRows (),
78  std::invalid_argument,
79  "Ifpack2::ReorderFilter: The input matrix is not square.");
80 
81  // Temp arrays for apply
82  Kokkos::resize(Indices_,A_->getLocalMaxNumRowEntries ());
83  Kokkos::resize(Values_,A_->getLocalMaxNumRowEntries ());
84 }
85 
86 
87 template<class MatrixType>
89 
90 
91 template<class MatrixType>
93 {
94  return A_->getComm();
95 }
96 
97 
98 
99 
100 template<class MatrixType>
103 {
105  A_.is_null (), std::runtime_error, "Ifpack2::ReorderFilter::"
106  "getRowMap: The matrix A is null, so there is no row Map.");
107 
108  return A_->getRowMap ();
109 }
110 
111 
112 template<class MatrixType>
115 {
117  A_.is_null (), std::runtime_error, "Ifpack2::ReorderFilter::"
118  "getColMap: The matrix A is null, so there is no column Map.");
119 
120  return A_->getColMap();
121 }
122 
123 
124 template<class MatrixType>
127 {
129  A_.is_null (), std::runtime_error, "Ifpack2::ReorderFilter::"
130  "getDomainMap: The matrix A is null, so there is no domain Map.");
131 
132  return A_->getDomainMap();
133 }
134 
135 
136 template<class MatrixType>
139 {
141  A_.is_null (), std::runtime_error, "Ifpack2::ReorderFilter::"
142  "getRangeMap: The matrix A is null, so there is no range Map.");
143 
144  return A_->getRangeMap();
145 }
146 
147 
148 template<class MatrixType>
149 Teuchos::RCP<const Tpetra::RowGraph<typename MatrixType::local_ordinal_type,
150  typename MatrixType::global_ordinal_type,
151  typename MatrixType::node_type> >
153 {
154  throw std::runtime_error("Ifpack2::ReorderFilter: does not support getGraph.");
155 }
156 
157 
158 template<class MatrixType>
160 {
161  return A_->getGlobalNumRows();
162 }
163 
164 
165 template<class MatrixType>
167 {
168  return A_->getGlobalNumCols();
169 }
170 
171 
172 template<class MatrixType>
174 {
175  return A_->getLocalNumRows();
176 }
177 
178 
179 template<class MatrixType>
181 {
182  return A_->getLocalNumCols();
183 }
184 
185 
186 template<class MatrixType>
187 typename MatrixType::global_ordinal_type ReorderFilter<MatrixType>::getIndexBase() const
188 {
189  return A_->getIndexBase();
190 }
191 
192 
193 template<class MatrixType>
195 {
196  return A_->getGlobalNumEntries();
197 }
198 
199 
200 template<class MatrixType>
202 {
203  return A_->getLocalNumEntries();
204 }
205 
206 template<class MatrixType>
207 typename MatrixType::local_ordinal_type ReorderFilter<MatrixType>::getBlockSize() const
208 {
209  return A_->getBlockSize();
210 }
211 
212 template<class MatrixType>
214 getNumEntriesInGlobalRow (global_ordinal_type globalRow) const
215 {
216  if (A_.is_null () || A_->getRowMap ().is_null ()) {
218  }
219  else {
220  const local_ordinal_type lclRow =
221  A_->getRowMap ()->getLocalElement (globalRow);
223  // The calling process doesn't own any entries in this row.
224  return static_cast<size_t> (0);
225  } else {
226  const local_ordinal_type origLclRow = reverseperm_[lclRow];
227  return A_->getNumEntriesInLocalRow (origLclRow);
228  }
229  }
230 }
231 
232 template<class MatrixType>
234 getNumEntriesInLocalRow (local_ordinal_type localRow) const
235 {
236  // Make sure that localRow is in bounds before using it to index
237  // into the permutation.
238  if (A_->getRowMap ()->isNodeLocalElement (localRow)) {
239  // localRow is a valid index into reverseperm_.
240  const local_ordinal_type localReorderedRow = reverseperm_[localRow];
241  return A_->getNumEntriesInLocalRow (localReorderedRow);
242  } else {
243  // The calling process doesn't own any entries in this row.
244  return static_cast<size_t> (0);
245  }
246 }
247 
248 
249 template<class MatrixType>
251 {
252  return A_->getGlobalMaxNumRowEntries();
253 }
254 
255 
256 template<class MatrixType>
258 {
259  return A_->getLocalMaxNumRowEntries();
260 }
261 
262 
263 template<class MatrixType>
265 {
266  return true;
267 }
268 
269 
270 template<class MatrixType>
272 {
273  return A_->isLocallyIndexed();
274 }
275 
276 
277 template<class MatrixType>
279 {
280  return A_->isGloballyIndexed();
281 }
282 
283 
284 template<class MatrixType>
286 {
287  return A_->isFillComplete();
288 }
289 
290 
291 template<class MatrixType>
293  getGlobalRowCopy (global_ordinal_type globalRow,
294  nonconst_global_inds_host_view_type &globalInd,
295  nonconst_values_host_view_type &val,
296  size_t& numEntries) const
297 {
298  using Teuchos::Array;
299  using Teuchos::ArrayView;
300  using Teuchos::av_reinterpret_cast;
301  typedef local_ordinal_type LO;
302  typedef Teuchos::OrdinalTraits<LO> OTLO;
303 
304  const map_type& rowMap = * (A_->getRowMap ());
305  const local_ordinal_type localRow = rowMap.getLocalElement (globalRow);
307  localRow == OTLO::invalid (), std::invalid_argument, "Ifpack2::Reorder"
308  "Filter::getGlobalRowCopy: The given global row index " << globalRow
309  << " is not owned by the calling process with rank "
310  << rowMap.getComm ()->getRank () << ".");
311 
312  // The Indices_ temp array is only used in apply, not getLocalRowCopy, so this is safe
313  numEntries = this->getNumEntriesInLocalRow (localRow);
314  this->getLocalRowCopy (localRow, Indices_, val, numEntries);
315 
316  // Convert local indices back to global indices.
317  for (size_t k = 0; k < numEntries; ++k) {
318  globalInd[k] = rowMap.getGlobalElement (Indices_[k]);
319  }
320 }
321 
322 
323 template<class MatrixType>
325 getLocalRowCopy (local_ordinal_type LocalRow,
326  nonconst_local_inds_host_view_type &Indices,
327  nonconst_values_host_view_type &Values,
328  size_t& NumEntries) const
329 
330 {
332  ! A_->getRowMap ()->isNodeLocalElement (LocalRow),
333  std::invalid_argument,
334  "Ifpack2::ReorderFilter::getLocalRowCopy: The given local row index "
335  << LocalRow << " is not a valid local row index on the calling process "
336  "with rank " << A_->getRowMap ()->getComm ()->getRank () << ".");
337 
338  // This duplicates code in getNumEntriesInGlobalRow, but avoids an
339  // extra array lookup and some extra tests.
340  const local_ordinal_type origLclRow = reverseperm_[LocalRow];
341  const size_t numEntries = A_->getNumEntriesInLocalRow (origLclRow);
342 
344  static_cast<size_t> (Indices.size ()) < numEntries ||
345  static_cast<size_t> (Values.size ()) < numEntries,
346  std::invalid_argument,
347  "Ifpack2::ReorderFilter::getLocalRowCopy: The given array views are not "
348  "long enough to store all the data in the given row " << LocalRow
349  << ". Indices.size() = " << Indices.size () << ", Values.size() = "
350  << Values.size () << ", but the (original) row has " << numEntries
351  << " entry/ies.");
352 
353  A_->getLocalRowCopy (origLclRow, Indices, Values, NumEntries);
354  // Do a col reindex via perm
355  //
356  // FIXME (mfh 30 Jan 2014) This assumes that the row and column
357  // indices are the same.
358  for (size_t i = 0; i < NumEntries; ++i) {
359  Indices[i] = perm_[Indices[i]];
360  }
361 }
362 
363 
364 template<class MatrixType>
365 void ReorderFilter<MatrixType>::getGlobalRowView(global_ordinal_type /* GlobalRow */,
366  global_inds_host_view_type &/*indices*/,
367  values_host_view_type &/*values*/) const
368 {
369  throw std::runtime_error("Ifpack2::ReorderFilter: does not support getGlobalRowView.");
370 }
371 
372 
373 
374 template<class MatrixType>
375 void ReorderFilter<MatrixType>::getLocalRowView(local_ordinal_type /* LocalRow */,
376  local_inds_host_view_type & /*indices*/,
377  values_host_view_type & /*values*/) const
378 {
379  throw std::runtime_error("Ifpack2::ReorderFilter: does not support getLocalRowView.");
380 }
381 
382 
383 
384 template<class MatrixType>
386 getLocalDiagCopy (Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &diag) const
387 {
388  // This is somewhat dubious as to how the maps match.
389  return A_->getLocalDiagCopy(diag);
390 }
391 
392 
393 template<class MatrixType>
394 void ReorderFilter<MatrixType>::leftScale(const Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& /* x */)
395 {
396  throw std::runtime_error("Ifpack2::ReorderFilter does not support leftScale.");
397 }
398 
399 
400 template<class MatrixType>
401 void ReorderFilter<MatrixType>::rightScale(const Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& /* x */)
402 {
403  throw std::runtime_error("Ifpack2::ReorderFilter does not support rightScale.");
404 }
405 
406 
407 template<class MatrixType>
409 apply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &X,
410  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &Y,
411  Teuchos::ETransp mode,
412  scalar_type alpha,
413  scalar_type beta) const
414 {
416 
418  alpha != STS::one () || beta != STS::zero (), std::logic_error,
419  "Ifpack2::ReorderFilter::apply is only implemented for alpha = 1 and "
420  "beta = 0. You set alpha = " << alpha << " and beta = " << beta << ".");
421 
422  // Note: This isn't AztecOO compliant. But neither was Ifpack's version.
423  // Note: The localized maps mean the matvec is trivial (and has no import)
425  X.getNumVectors() != Y.getNumVectors(), std::runtime_error,
426  "Ifpack2::ReorderFilter::apply: X.getNumVectors() != Y.getNumVectors().");
427 
428  const scalar_type zero = STS::zero ();
430  Teuchos::ArrayRCP<Teuchos::ArrayRCP<scalar_type> > y_ptr = Y.get2dViewNonConst();
431 
432  Y.putScalar (zero);
433  const size_t NumVectors = Y.getNumVectors ();
434 
435  for (size_t i = 0; i < A_->getLocalNumRows (); ++i) {
436  size_t Nnz;
437  // Use this class's getrow to make the below code simpler
438  getLocalRowCopy (i, Indices_ , Values_ , Nnz);
439  scalar_type* Values = reinterpret_cast<scalar_type*>(Values_.data());
440  if (mode == Teuchos::NO_TRANS) {
441  for (size_t j = 0; j < Nnz; ++j) {
442  for (size_t k = 0; k < NumVectors; ++k) {
443  y_ptr[k][i] += Values[j] * x_ptr[k][Indices_[j]];
444  }
445  }
446  }
447  else if (mode == Teuchos::TRANS) {
448  for (size_t j = 0; j < Nnz; ++j) {
449  for (size_t k = 0; k < NumVectors; ++k) {
450  y_ptr[k][Indices_[j]] += Values[j] * x_ptr[k][i];
451  }
452  }
453  }
454  else { //mode==Teuchos::CONJ_TRANS
455  for (size_t j = 0; j < Nnz; ++j) {
456  for (size_t k = 0; k < NumVectors; ++k) {
457  y_ptr[k][Indices_[j]] += STS::conjugate(Values[j]) * x_ptr[k][i];
458  }
459  }
460  }
461  }
462 }
463 
464 
465 template<class MatrixType>
467 {
468  return true;
469 }
470 
471 
472 template<class MatrixType>
474 {
475  return false;
476 }
477 
478 
479 template<class MatrixType>
480 typename ReorderFilter<MatrixType>::mag_type ReorderFilter<MatrixType>::getFrobeniusNorm() const
481 {
482  // Reordering doesn't change the Frobenius norm.
483  return A_->getFrobeniusNorm ();
484 }
485 
486 
487 template<class MatrixType>
489 permuteOriginalToReordered (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &originalX,
490  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &reorderedY) const
491 {
492  this->template permuteOriginalToReorderedTempl<scalar_type,scalar_type>(originalX, reorderedY);
493 }
494 
495 
496 template<class MatrixType>
497 template<class DomainScalar, class RangeScalar>
498 void ReorderFilter<MatrixType>::permuteOriginalToReorderedTempl(const Tpetra::MultiVector<DomainScalar,local_ordinal_type,global_ordinal_type,node_type> &originalX,
499  Tpetra::MultiVector<RangeScalar,local_ordinal_type,global_ordinal_type,node_type> &reorderedY) const
500 {
501  TEUCHOS_TEST_FOR_EXCEPTION(originalX.getNumVectors() != reorderedY.getNumVectors(), std::runtime_error,
502  "Ifpack2::ReorderFilter::permuteOriginalToReordered ERROR: X.getNumVectors() != Y.getNumVectors().");
503 
504  Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > x_ptr = originalX.get2dView();
505  Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > y_ptr = reorderedY.get2dViewNonConst();
506 
507  const local_ordinal_type blockSize = getBlockSize();
508  const local_ordinal_type numRows = originalX.getLocalLength() / blockSize;
509  for(size_t k=0; k < originalX.getNumVectors(); k++)
510  for(local_ordinal_type i=0; i< numRows; i++)
511  for(local_ordinal_type j=0; j< blockSize; ++j)
512  y_ptr[k][perm_[i]*blockSize + j] = (RangeScalar)x_ptr[k][i*blockSize + j];
513 }
514 
515 
516 template<class MatrixType>
517 void ReorderFilter<MatrixType>::permuteReorderedToOriginal(const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &reorderedX,
518  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &originalY) const
519 {
520  this->template permuteReorderedToOriginalTempl<scalar_type,scalar_type>(reorderedX, originalY);
521 }
522 
523 
524 template<class MatrixType>
525 template<class DomainScalar, class RangeScalar>
527 permuteReorderedToOriginalTempl (const Tpetra::MultiVector<DomainScalar,local_ordinal_type,global_ordinal_type,node_type> &reorderedX,
528  Tpetra::MultiVector<RangeScalar,local_ordinal_type,global_ordinal_type,node_type> &originalY) const
529 {
531  reorderedX.getNumVectors() != originalY.getNumVectors(),
532  std::runtime_error,
533  "Ifpack2::ReorderFilter::permuteReorderedToOriginal: "
534  "X.getNumVectors() != Y.getNumVectors().");
535 
536 #ifdef HAVE_IFPACK2_DEBUG
537  {
539  Teuchos::Array<magnitude_type> norms (reorderedX.getNumVectors ());
540  reorderedX.norm2 (norms ());
541  bool good = true;
542  for (size_t j = 0;
543  j < reorderedX.getNumVectors (); ++j) {
544  if (STM::isnaninf (norms[j])) {
545  good = false;
546  break;
547  }
548  }
550  ! good, std::runtime_error, "Ifpack2::ReorderFilter::"
551  "permuteReorderedToOriginalTempl: The 2-norm of the input reorderedX is "
552  "NaN or Inf.");
553  }
554 #endif // HAVE_IFPACK2_DEBUG
555 
556  Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > x_ptr = reorderedX.get2dView();
557  Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > y_ptr = originalY.get2dViewNonConst();
558 
559  const local_ordinal_type blockSize = getBlockSize();
560  const local_ordinal_type numRows = reorderedX.getLocalLength() / blockSize;
561  for (size_t k = 0; k < reorderedX.getNumVectors (); ++k) {
562  for (local_ordinal_type i = 0; i < numRows; ++i) {
563  for(local_ordinal_type j = 0; j < blockSize; ++j) {
564  y_ptr[k][reverseperm_[i]*blockSize + j] = (RangeScalar) x_ptr[k][i*blockSize + j];
565  }
566  }
567  }
568 
569 #ifdef HAVE_IFPACK2_DEBUG
570  {
572  Teuchos::Array<magnitude_type> norms (originalY.getNumVectors ());
573  originalY.norm2 (norms ());
574  bool good = true;
575  for (size_t j = 0;
576  j < originalY.getNumVectors (); ++j) {
577  if (STM::isnaninf (norms[j])) {
578  good = false;
579  break;
580  }
581  }
583  ! good, std::runtime_error, "Ifpack2::ReorderFilter::"
584  "permuteReorderedToOriginalTempl: The 2-norm of the output originalY is "
585  "NaN or Inf.");
586  }
587 #endif // HAVE_IFPACK2_DEBUG
588 }
589 
590 } // namespace Ifpack2
591 
592 #define IFPACK2_REORDERFILTER_INSTANT(S,LO,GO,N) \
593  template class Ifpack2::ReorderFilter< Tpetra::RowMatrix<S, LO, GO, N> >;
594 
595 #endif
virtual size_t getLocalNumRows() const
Returns the number of rows owned on the calling node.
Definition: Ifpack2_ReorderFilter_def.hpp:173
Wraps a Tpetra::RowMatrix in a filter that reorders local rows and columns.
Definition: Ifpack2_ReorderFilter_decl.hpp:69
virtual mag_type getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:480
virtual Teuchos::RCP< const map_type > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:102
virtual bool supportsRowViews() const
Returns true if RowViews are supported.
Definition: Ifpack2_ReorderFilter_def.hpp:473
virtual size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
Definition: Ifpack2_ReorderFilter_def.hpp:250
virtual Teuchos::RCP< const map_type > getDomainMap() const
Returns the Map that describes the domain distribution in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:126
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:159
virtual size_t getLocalNumCols() const
Returns the number of columns needed to apply the forward operator on this node, i.e., the number of elements listed in the column map.
Definition: Ifpack2_ReorderFilter_def.hpp:180
virtual Teuchos::RCP< const map_type > getColMap() const
Returns the Map that describes the column distribution in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:114
ReorderFilter(const Teuchos::RCP< const row_matrix_type > &A, const Teuchos::ArrayRCP< local_ordinal_type > &perm, const Teuchos::ArrayRCP< local_ordinal_type > &reverseperm)
Constructor.
Definition: Ifpack2_ReorderFilter_def.hpp:58
virtual Teuchos::RCP< const map_type > getRangeMap() const
Returns the Map that describes the range distribution in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:138
virtual void getGlobalRowView(global_ordinal_type 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_ReorderFilter_def.hpp:365
virtual void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
, where Op(A) is either A, , or .
Definition: Ifpack2_ReorderFilter_def.hpp:409
virtual void getGlobalRowCopy(global_ordinal_type 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_ReorderFilter_def.hpp:293
virtual bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
Definition: Ifpack2_ReorderFilter_def.hpp:278
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The matrix&#39;s communicator.
Definition: Ifpack2_ReorderFilter_def.hpp:92
virtual void getLocalDiagCopy(Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
Definition: Ifpack2_ReorderFilter_def.hpp:386
virtual size_t getNumEntriesInLocalRow(local_ordinal_type localRow) const
The current number of entries in this matrix, stored on the calling process, in the row whose local i...
Definition: Ifpack2_ReorderFilter_def.hpp:234
virtual void getLocalRowView(local_ordinal_type 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_ReorderFilter_def.hpp:375
virtual bool hasColMap() const
Indicates whether this matrix has a well-defined column map.
Definition: Ifpack2_ReorderFilter_def.hpp:264
virtual size_t getNumEntriesInGlobalRow(global_ordinal_type globalRow) const
The current number of entries in this matrix, stored on the calling process, in the row whose global ...
Definition: Ifpack2_ReorderFilter_def.hpp:214
virtual global_size_t getGlobalNumCols() const
Returns the number of global columns in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:166
virtual bool hasTransposeApply() const
Whether apply() can apply the transpose or conjugate transpose.
Definition: Ifpack2_ReorderFilter_def.hpp:466
virtual size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:201
virtual global_ordinal_type getIndexBase() const
Returns the index base for global indices for this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:187
virtual global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:194
virtual size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
Definition: Ifpack2_ReorderFilter_def.hpp:257
virtual bool isFillComplete() const
Returns true if fillComplete() has been called.
Definition: Ifpack2_ReorderFilter_def.hpp:285
virtual void permuteReorderedToOriginal(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &reorderedX, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &originalY) const
Permute multivector: reordered-to-original.
Definition: Ifpack2_ReorderFilter_def.hpp:517
virtual Teuchos::RCP< const Tpetra::RowGraph< local_ordinal_type, global_ordinal_type, node_type > > getGraph() const
Returns the RowGraph associated with this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:152
virtual void leftScale(const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &x)
Scales the RowMatrix on the left with the Vector x.
Definition: Ifpack2_ReorderFilter_def.hpp:394
virtual void rightScale(const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &x)
Scales the RowMatrix on the right with the Vector x.
Definition: Ifpack2_ReorderFilter_def.hpp:401
virtual void getLocalRowCopy(local_ordinal_type 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_ReorderFilter_def.hpp:325
virtual void permuteOriginalToReordered(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &originalX, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &reorderedY) const
Permute multivector: original-to-reordered.
Definition: Ifpack2_ReorderFilter_def.hpp:489
virtual ~ReorderFilter()
Destructor.
Definition: Ifpack2_ReorderFilter_def.hpp:88
virtual bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
Definition: Ifpack2_ReorderFilter_def.hpp:271
virtual local_ordinal_type getBlockSize() const
The number of degrees of freedom per mesh point.
Definition: Ifpack2_ReorderFilter_def.hpp:207
bool is_null() const