Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_CrsMatrixWrap_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Xpetra: A linear algebra interface package
4//
5// Copyright 2012 NTESS and the Xpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10// WARNING: This code is experimental. Backwards compatibility should not be expected.
11
12#ifndef XPETRA_CRSMATRIXWRAP_DECL_HPP
13#define XPETRA_CRSMATRIXWRAP_DECL_HPP
14
15#include <Tpetra_KokkosCompat_DefaultNode.hpp>
16
17#include "Xpetra_ConfigDefs.hpp"
18#include "Xpetra_Exceptions.hpp"
19
21#include "Xpetra_CrsGraph.hpp"
22#include "Xpetra_CrsMatrix.hpp"
24
25#include "Xpetra_Matrix.hpp"
27#include "Xpetra_TpetraBlockCrsMatrix.hpp"
29#include "Xpetra_TpetraCrsMatrix.hpp"
30
32#include <Teuchos_Hashtable.hpp>
33
34#include <Tpetra_Filter.hpp>
35
40namespace Xpetra {
41
42typedef std::string viewLabel_t;
43
48template <class Scalar,
49 class LocalOrdinal,
50 class GlobalOrdinal,
51 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
52class CrsMatrixWrap : public Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
63
64 public:
66
67
69 CrsMatrixWrap(const RCP<const Map> &rowMap);
70
72 CrsMatrixWrap(const RCP<const Map> &rowMap,
73 size_t maxNumEntriesPerRow);
74
76 CrsMatrixWrap(const RCP<const Map> &rowMap,
77 const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
78
80 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, size_t maxNumEntriesPerRow);
81
83 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
84
86 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const local_matrix_type &lclMatrix, const Teuchos::RCP<Teuchos::ParameterList> &params = null);
87
89 CrsMatrixWrap(const local_matrix_type &lclMatrix, const RCP<const Map> &rowMap, const RCP<const Map> &colMap,
90 const RCP<const Map> &domainMap = Teuchos::null, const RCP<const Map> &rangeMap = Teuchos::null,
91 const Teuchos::RCP<Teuchos::ParameterList> &params = null);
92
94
95 CrsMatrixWrap(const RCP<const CrsGraph> &graph, const RCP<ParameterList> &paramList = Teuchos::null);
96
99 const RCP<ParameterList> &paramList = Teuchos::null);
100
102 virtual ~CrsMatrixWrap();
103
105
107
108
110
121 void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals);
122
124
131 void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals);
132
134
139 void replaceGlobalValues(GlobalOrdinal globalRow,
141 const ArrayView<const Scalar> &vals);
142
144
147 void replaceLocalValues(LocalOrdinal localRow,
149 const ArrayView<const Scalar> &vals);
150
152 virtual void setAllToScalar(const Scalar &alpha);
153
155 void scale(const Scalar &alpha);
156
158
160
161
170 void resumeFill(const RCP<ParameterList> &params = null);
171
183 void fillComplete(const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap, const RCP<Teuchos::ParameterList> &params = null);
184
198 // TODO : Get ride of "Tpetra"::OptimizeOption
199 void fillComplete(const RCP<ParameterList> &params = null);
200
202
204
207
209
212
214 size_t getLocalNumRows() const;
215
218
220 size_t getLocalNumEntries() const;
221
223
224 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
225
227
228 size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
229
231
233 size_t getGlobalMaxNumRowEntries() const;
234
236
238 size_t getLocalMaxNumRowEntries() const;
239
241 bool isLocallyIndexed() const;
242
244 bool isGloballyIndexed() const;
245
247 bool isFillComplete() const;
248
250
262 void getLocalRowCopy(LocalOrdinal LocalRow,
263 const ArrayView<LocalOrdinal> &Indices,
264 const ArrayView<Scalar> &Values,
265 size_t &NumEntries) const;
266
268
277 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) const;
278
280
289 void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) const;
290
292
295
298
301
304
307
310
312 bool haveGlobalConstants() const;
313
315
317
318
320
329 // TODO virtual=0 // TODO: Add default parameters ?
330 // void multiply(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & X, MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y, Teuchos::ETransp trans, Scalar alpha, Scalar beta) const {
331 // matrixData_->multiply(X, Y, trans, alpha, beta);
332 // }
333
335
337
338
340
346 Scalar alpha = ScalarTraits<Scalar>::one(),
347 Scalar beta = ScalarTraits<Scalar>::zero()) const;
348
352 Teuchos::ETransp mode,
353 Scalar alpha,
354 Scalar beta,
355 bool sumInterfaceValues,
356 const RCP<Import<LocalOrdinal, GlobalOrdinal, Node>> &regionInterfaceImporter,
357 const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const;
358
362
366
369 const RCP<const Map> &getColMap() const;
370
372 const RCP<const Map> &getColMap(viewLabel_t viewLabel) const;
373
375
377
379 //{@
380
383
385 void doImport(const Matrix &source,
387
389 void doExport(const Matrix &dest,
391
393 void doImport(const Matrix &source,
395
397 void doExport(const Matrix &dest,
399
400 // @}
401
403
404
406 std::string description() const;
407
410
412
413 void setObjectLabel(const std::string &objectLabel);
415
417 virtual typename local_matrix_type::host_mirror_type getLocalMatrixHost() const;
418
419 // JG: Added:
420
421 bool hasCrsGraph() const;
422
425
427
429 LocalOrdinal GetStorageBlockSize() const;
430
435
438
440 private:
441 // Default view is created after fillComplete()
442 // Because ColMap might not be available before fillComplete().
443 void CreateDefaultView();
444
445 // The colMap can be <tt>null</tt> until fillComplete() is called. The default view of the Matrix have to be updated when fillComplete() is called.
446 // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated when getColMap() is called.
447 void updateDefaultView() const;
448 // The boolean finalDefaultView_ keep track of the status of the default view (= already updated or not)
449 // See also CrsMatrixWrap::updateDefaultView()
450 mutable bool finalDefaultView_;
451
452 // The underlying matrix object
454
455}; // class CrsMatrixWrap
456
457template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
460 return Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(matrix, true)->getCrsMatrix();
461}
462
463template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
466 return Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(matrix, true)->getCrsMatrix();
467}
468
469template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
474
475template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
480
481template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
485 auto xpCrs = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
487}
488
489template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
493 auto xpCrs = Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
495}
496
497template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
501 auto xpCrs = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
503}
504
505template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
509 auto xpCrs = Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
511}
512
513template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
518
519template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
524
525template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
530
531template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
536
537template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
540 return Teuchos::rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true)->getTpetra_BlockCrsMatrix();
541}
542
543template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
546 return Teuchos::rcp_dynamic_cast<Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true)->getTpetra_BlockCrsMatrixNonConst();
547}
548
549template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
554
555template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
560
561template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
566
567template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
572
573template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
578
579template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
584
585template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
588 auto crsMat = toCrsMatrix(mat);
589 auto tmp_Crs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
590 if (!tmp_Crs.is_null()) {
591 return tmp_Crs->getTpetra_CrsMatrixNonConst();
592 } else {
593 auto tmp_BlockCrs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
594 if (tmp_BlockCrs.is_null())
595 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
596 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
597 }
598}
599
600template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
603 auto crsMat = toCrsMatrix(mat);
604 auto tmp_Crs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
605 if (!tmp_Crs.is_null()) {
606 return tmp_Crs->getTpetra_CrsMatrix();
607 } else {
608 auto tmp_BlockCrs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
609 if (tmp_BlockCrs.is_null())
610 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
611 return tmp_BlockCrs->getTpetra_BlockCrsMatrix();
612 }
613}
614
615template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
618 auto mat = Teuchos::rcp_dynamic_cast<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
619 auto rmat = Teuchos::rcp_dynamic_cast<Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
620 if (!mat.is_null()) {
621 return toTpetraRowMatrix(mat);
622 } else if (!rmat.is_null()) {
623 return rmat->getTpetra_RowMatrixNonConst();
624 } else {
625 auto tpOp = Teuchos::rcp_dynamic_cast<Xpetra::TpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true);
626 auto tOp = tpOp->getOperator();
627 auto tRow = Teuchos::rcp_dynamic_cast<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp, true);
628 return tRow;
629 }
630}
631
632template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
635 auto mat = Teuchos::rcp_dynamic_cast<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
636 auto rmat = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
637 if (!mat.is_null()) {
638 return toTpetraRowMatrix(mat);
639 } else if (!rmat.is_null()) {
640 return rmat->getTpetra_RowMatrix();
641 } else {
642 auto tpOp = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true);
643 auto tOp = tpOp->getOperator();
644 auto tRow = Teuchos::rcp_dynamic_cast<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp, true);
645 return tRow;
646 }
647}
648
649template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
652 auto tpCrsA = Teuchos::rcp_dynamic_cast<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsA);
653 if (!tpCrsA.is_null()) {
654 auto tpCrsFiltered = Tpetra::applyFilter_vals(*tpCrsA, filter);
655 crsFiltered = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tpCrsFiltered);
656 }
657 return crsFiltered;
658}
659
660template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
662 auto crsA = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A)->getCrsMatrix();
663 auto crsFiltered = applyFilter_vals(crsA, filter);
664 auto filtered = rcp(new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(crsFiltered));
665 return filtered;
666}
667
668template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
671 auto tpCrsA = Teuchos::rcp_dynamic_cast<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsA);
672 if (!tpCrsA.is_null()) {
673 auto tpCrsFiltered = Tpetra::applyFilter_LID(*tpCrsA, filter);
674 crsFiltered = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tpCrsFiltered);
675 }
676 return crsFiltered;
677}
678
679template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
681 auto crsA = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A)->getCrsMatrix();
682 auto crsFiltered = applyFilter_LID(crsA, filter);
683 auto filtered = rcp(new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(crsFiltered));
684 return filtered;
685}
686
687template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
690 auto tpCrsA = Teuchos::rcp_dynamic_cast<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsA);
691 if (!tpCrsA.is_null()) {
692 auto tpCrsFiltered = Tpetra::applyFilter_GID(*tpCrsA, filter);
693 crsFiltered = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tpCrsFiltered);
694 }
695 return crsFiltered;
696}
697
698template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
700 auto crsA = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A)->getCrsMatrix();
701 auto crsFiltered = applyFilter_GID(crsA, filter);
702 auto filtered = rcp(new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(crsFiltered));
703 return filtered;
704}
705
706} // namespace Xpetra
707
708#define XPETRA_CRSMATRIXWRAP_SHORT
709#endif // XPETRA_CRSMATRIXWRAP_DECL_HPP
710
711// NOTE: if CrsMatrix and VbrMatrix share a common interface for fillComplete() etc, I can move some stuff in Xpetra_Matrix.hpp
712// TODO: getUnderlyingMatrix() method
static const EVerbosityLevel verbLevel_default
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, execution_space, void, typename local_graph_type::size_type > local_matrix_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, execution_space, void, typename local_graph_type::size_type > local_matrix_device_type
typename local_matrix_device_type::host_mirror_type local_matrix_host_type
Concrete implementation of Xpetra::Matrix.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale matrix using the given vector entries.
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_host_type local_matrix_host_type
LocalOrdinal GetStorageBlockSize() const
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
void residual(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const
Compute a residual R = B - (*this) * X.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the specified global row.
RCP< const CrsGraph > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Matrix
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< Teuchos::ParameterList > &params=null)
Signal that data entry is complete, specifying domain and range maps.
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale matrix using the given vector entries.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Get Frobenius norm of the matrix.
std::string description() const
Return a simple one-line description of this object.
Xpetra::CrsMatrixFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixFactory
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
Xpetra::TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrix
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_device_type local_matrix_device_type
size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
bool hasCrsGraph() const
Supports the getCrsGraph() call.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type local_matrix_type
void resumeFill(const RCP< ParameterList > &params=null)
global_size_t getGlobalNumCols() const
Returns the number of global columns in the matrix.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
void setObjectLabel(const std::string &objectLabel)
virtual local_matrix_type::host_mirror_type getLocalMatrixHost() const
const RCP< const Map > & getColMap() const
Returns the Map that describes the column distribution in this matrix. This might be null until fillC...
const RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
Xpetra::MatrixView< Scalar, LocalOrdinal, GlobalOrdinal, Node > MatrixView
virtual local_matrix_type getLocalMatrixDevice() const
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void replaceCrsMatrix(RCP< CrsMatrix > &M)
Expert only.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void doExport(const Matrix &dest, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void getLocalDiagCopy(Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row idices.
virtual ~CrsMatrixWrap()
Destructor.
const RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
size_t getLocalNumRows() const
Returns the number of matrix rows owned on the calling node.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
void doImport(const Matrix &source, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
RCP< CrsMatrix > getCrsMatrix() const
virtual void apply(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode.
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > CrsGraph
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
Exception indicating invalid cast attempted.
Xpetra-specific matrix class.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< matrix_type > applyFilter_GID(const matrix_type &A, const filter_type &filter)
Teuchos::RCP< matrix_type > applyFilter_LID(const matrix_type &A, const filter_type &filter)
Teuchos::RCP< matrix_type > applyFilter_vals(const matrix_type &A, const filter_type &filter)
Teuchos::RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > toCrsMatrix(const Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &matrix)
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > graph)
size_t global_size_t
Global size_t object.
RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > applyFilter_vals(const RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &crsA, const filter_type &filter)
Teuchos::RCP< Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > toTpetraBlock(const Teuchos::RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A)
RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > applyFilter_LID(const RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &crsA, const filter_type &filter)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > applyFilter_GID(const RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &crsA, const filter_type &filter)
Teuchos::RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > toTpetraRowMatrix(const Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &mat)
CombineMode
Xpetra::Combine Mode enumerable type.