Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Related Functions | List of all members
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Class Template Referenceabstract

Sparse matrix that presents a row-oriented interface that lets users read or modify entries. More...

#include <Tpetra_CrsMatrix_decl.hpp>

Inheritance diagram for Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >:
Inheritance graph
[legend]

Public Types

Typedefs
using scalar_type = Scalar
 The type of each entry in the matrix. More...
 
using local_ordinal_type = LocalOrdinal
 The type of each local index in the matrix. More...
 
using global_ordinal_type = GlobalOrdinal
 The type of each global index in the matrix. More...
 
using device_type = typename Node::device_type
 The Kokkos device type. More...
 
using execution_space = typename device_type::execution_space
 The Kokkos execution space. More...
 
using memory_space = typename device_type::memory_space
 The Kokkos memory space. More...
 
using node_type = Node
 This class' Kokkos Node type. More...
 
using map_type = Map< LocalOrdinal, GlobalOrdinal, Node >
 The Map specialization suitable for this CrsMatrix specialization. More...
 
using import_type = Import< LocalOrdinal, GlobalOrdinal, Node >
 The Import specialization suitable for this CrsMatrix specialization. More...
 
using export_type = Export< LocalOrdinal, GlobalOrdinal, Node >
 The Export specialization suitable for this CrsMatrix specialization. More...
 
using row_matrix_type = RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >
 The RowMatrix representing the base class of CrsMatrix. More...
 
using impl_scalar_type = typename row_matrix_type::impl_scalar_type
 The type used internally in place of Scalar. More...
 
using mag_type = typename Kokkos::ArithTraits< impl_scalar_type >::mag_type
 Type of a norm result. More...
 
using crs_graph_type = CrsGraph< LocalOrdinal, GlobalOrdinal, Node >
 The CrsGraph specialization suitable for this CrsMatrix specialization. More...
 
using local_graph_device_type = typename crs_graph_type::local_graph_device_type
 The part of the sparse matrix's graph on each MPI process. More...
 
using local_graph_host_type = typename crs_graph_type::local_graph_host_type
 
using local_matrix_device_type = KokkosSparse::CrsMatrix< impl_scalar_type, local_ordinal_type, device_type, void, typename local_graph_device_type::size_type >
 The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI process. More...
 
using local_matrix_host_type = typename local_matrix_device_type::HostMirror
 
using local_multiply_op_type = LocalCrsMatrixOperator< scalar_type, scalar_type, device_type >
 The type of the local matrix-vector operator (a wrapper of KokkosSparse::CrsMatrix ) More...
 
using row_ptrs_device_view_type = typename row_matrix_type::row_ptrs_device_view_type
 
using row_ptrs_host_view_type = typename row_matrix_type::row_ptrs_host_view_type
 
using local_inds_device_view_type = typename row_matrix_type::local_inds_device_view_type
 
using local_inds_host_view_type = typename row_matrix_type::local_inds_host_view_type
 
using nonconst_local_inds_host_view_type = typename row_matrix_type::nonconst_local_inds_host_view_type
 
using global_inds_device_view_type = typename row_matrix_type::global_inds_device_view_type
 
using global_inds_host_view_type = typename row_matrix_type::global_inds_host_view_type
 
using nonconst_global_inds_host_view_type = typename row_matrix_type::nonconst_global_inds_host_view_type
 
using values_device_view_type = typename row_matrix_type::values_device_view_type
 
using values_host_view_type = typename row_matrix_type::values_host_view_type
 
using nonconst_values_host_view_type = typename row_matrix_type::nonconst_values_host_view_type
 
Typedefs
using packet_type = typename::Kokkos::ArithTraits< char >::val_type
 The type of each datum being sent or received in an Import or Export. More...
 

Public Member Functions

local_matrix_host_type::values_type::const_type getLocalValuesHost (Access::ReadOnlyStruct s) const
 Get the Kokkos local values on host, read only. More...
 
local_matrix_host_type::values_type getLocalValuesHost (Access::ReadWriteStruct s)
 Get the Kokkos local values on host, read write. More...
 
local_matrix_host_type::values_type getLocalValuesHost (Access::OverwriteAllStruct s)
 Get the Kokkos local values on host, overwrite all. More...
 
local_matrix_device_type::values_type::const_type getLocalValuesDevice (Access::ReadOnlyStruct s) const
 Get the Kokkos local values on device, read only. More...
 
local_matrix_device_type::values_type getLocalValuesDevice (Access::ReadWriteStruct s)
 Get the Kokkos local values on device, read write. More...
 
local_matrix_device_type::values_type getLocalValuesDevice (Access::OverwriteAllStruct s)
 Get the Kokkos local values on device, overwrite all. More...
 
void importAndFillComplete (Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &destMatrix, const import_type &importer, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) const
 Import from this to the given destination matrix, and make the result fill complete. More...
 
void importAndFillComplete (Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &destMatrix, const import_type &rowImporter, const import_type &domainImporter, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params) const
 Import from this to the given destination matrix, and make the result fill complete. More...
 
void exportAndFillComplete (Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &destMatrix, const export_type &exporter, const Teuchos::RCP< const map_type > &domainMap=Teuchos::null, const Teuchos::RCP< const map_type > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) const
 Export from this to the given destination matrix, and make the result fill complete. More...
 
void exportAndFillComplete (Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &destMatrix, const export_type &rowExporter, const export_type &domainExporter, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params) const
 Export from this to the given destination matrix, and make the result fill complete. More...
 
bool haveGlobalConstants () const
 Returns true if globalConstants have been computed; false otherwise. More...
 
Transformational methods
void globalAssemble ()
 Communicate nonlocal contributions to other processes. More...
 
void resumeFill (const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Resume operations that may change the values or structure of the matrix. More...
 
void fillComplete (const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Tell the matrix that you are done changing its structure or values, and that you are ready to do computational kernels (e.g., sparse matrix-vector multiply) with it. More...
 
void fillComplete (const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Tell the matrix that you are done changing its structure or values, and that you are ready to do computational kernels (e.g., sparse matrix-vector multiply) with it. Set default domain and range Maps. More...
 
void expertStaticFillComplete (const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< const import_type > &importer=Teuchos::null, const Teuchos::RCP< const export_type > &exporter=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Perform a fillComplete on a matrix that already has data. More...
 
void replaceColMap (const Teuchos::RCP< const map_type > &newColMap)
 Replace the matrix's column Map with the given Map. More...
 
void reindexColumns (crs_graph_type *const graph, const Teuchos::RCP< const map_type > &newColMap, const Teuchos::RCP< const import_type > &newImport=Teuchos::null, const bool sortEachRow=true)
 Reindex the column indices in place, and replace the column Map. Optionally, replace the Import object as well. More...
 
void replaceDomainMap (const Teuchos::RCP< const map_type > &newDomainMap)
 Replace the current domain Map with the given objects. More...
 
void replaceDomainMapAndImporter (const Teuchos::RCP< const map_type > &newDomainMap, Teuchos::RCP< const import_type > &newImporter)
 Replace the current domain Map and Import with the given objects. More...
 
void replaceRangeMap (const Teuchos::RCP< const map_type > &newRangeMap)
 Replace the current range Map with the given objects. More...
 
void replaceRangeMapAndExporter (const Teuchos::RCP< const map_type > &newRangeMap, Teuchos::RCP< const export_type > &newExporter)
 Replace the current Range Map and Export with the given objects. More...
 
virtual void removeEmptyProcessesInPlace (const Teuchos::RCP< const map_type > &newMap) override
 Remove processes owning zero rows from the Maps and their communicator. More...
 
Local apply
void localApply (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, const Teuchos::ETransp mode=Teuchos::NO_TRANS, const Scalar &alpha=Teuchos::ScalarTraits< Scalar >::one(), const Scalar &beta=Teuchos::ScalarTraits< Scalar >::zero()) const
 Compute the local part of a sparse matrix-(Multi)Vector multiply. More...
 
template<class T >
Teuchos::RCP< CrsMatrix< T,
LocalOrdinal, GlobalOrdinal,
Node > > 
convert () const
 Return another CrsMatrix with the same entries, but converted to a different Scalar type T. More...
 
Methods implementing Operator
void apply (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, 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 override
 Compute a sparse matrix-MultiVector multiply. More...
 
bool hasTransposeApply () const override
 Whether apply() allows applying the transpose or conjugate transpose. More...
 
Teuchos::RCP< const map_typegetDomainMap () const override
 The domain Map of this matrix. More...
 
Teuchos::RCP< const map_typegetRangeMap () const override
 The range Map of this matrix. More...
 
Other "apply"-like methods
virtual Teuchos::RCP
< RowMatrix< Scalar,
LocalOrdinal, GlobalOrdinal,
Node > > 
add (const Scalar &alpha, const RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params) const override
 Implementation of RowMatrix::add: return alpha*A + beta*this. More...
 
Implementation of Teuchos::Describable interface
std::string description () const override
 A one-line description of this object. More...
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
 Print this object with the given verbosity level to the given output stream. More...
 
Extraction Methods
virtual void getGlobalRowCopy (GlobalOrdinal GlobalRow, nonconst_global_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const =0
 Get a copy of the given global row's entries. More...
 
virtual void getLocalRowCopy (LocalOrdinal LocalRow, nonconst_local_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const =0
 Get a copy of the given local row's entries. More...
 
virtual void getGlobalRowView (GlobalOrdinal GlobalRow, global_inds_host_view_type &indices, values_host_view_type &values) const =0
 Get a constant, nonpersisting, globally indexed view of the given row of the matrix. More...
 
virtual void getLocalRowView (LocalOrdinal LocalRow, local_inds_host_view_type &indices, values_host_view_type &values) const =0
 Get a constant, nonpersisting, locally indexed view of the given row of the matrix. More...
 
Implementation of Packable interface
virtual void pack (const Teuchos::ArrayView< const LocalOrdinal > &exportLIDs, Teuchos::Array< char > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets) const
 Pack this object's data for an Import or Export. More...
 
Pure virtual functions to be overridden by subclasses.
virtual bool hasDiagonal () const
 Whether this operator can return its diagonal. More...
 
Public methods for redistributing data
void doImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 Import data into this object using an Import object ("forward mode"). More...
 
void doImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 Import data into this object using an Export object ("reverse mode"). More...
 
void doExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 Export data into this object using an Export object ("forward mode"). More...
 
void doExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 Export data into this object using an Import object ("reverse mode"). More...
 
void beginImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 
void beginImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void beginExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void beginExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 
void endImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 
void endImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void endExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void endExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 
bool transferArrived () const
 Whether the data from an import/export operation has arrived, and is ready for the unpack and combine step. More...
 
Attribute accessor methods
bool isDistributed () const
 Whether this is a globally distributed object. More...
 
virtual Teuchos::RCP< const
map_type
getMap () const
 The Map describing the parallel distribution of this object. More...
 
I/O methods
void print (std::ostream &os) const
 Print this object to the given output stream. More...
 
Methods for use only by experts
virtual void removeEmptyProcessesInPlace (const Teuchos::RCP< const map_type > &newMap)
 Remove processes which contain no entries in this object's Map. More...
 

Protected Member Functions

virtual void insertGlobalValuesImpl (crs_graph_type &graph, RowInfo &rowInfo, const GlobalOrdinal gblColInds[], const impl_scalar_type vals[], const size_t numInputEnt)
 Common implementation detail of insertGlobalValues and insertGlobalValuesFiltered. More...
 
void allocateValues (ELocalGlobal lg, GraphAllocationStatus gas, const bool verbose)
 Allocate values (and optionally indices) using the Node. More...
 
void sortAndMergeIndicesAndValues (const bool sorted, const bool merged)
 Sort and merge duplicate local column indices in all rows on the calling process, along with their corresponding values. More...
 
Teuchos::RCP< MVgetColumnMapMultiVector (const MV &X_domainMap, const bool force=false) const
 Create a (or fetch a cached) column Map MultiVector. More...
 
Teuchos::RCP< MVgetRowMapMultiVector (const MV &Y_rangeMap, const bool force=false) const
 Create a (or fetch a cached) row Map MultiVector. More...
 
void applyNonTranspose (const MV &X_in, MV &Y_in, Scalar alpha, Scalar beta) const
 Special case of apply() for mode == Teuchos::NO_TRANS. More...
 
void applyTranspose (const MV &X_in, MV &Y_in, const Teuchos::ETransp mode, Scalar alpha, Scalar beta) const
 Special case of apply() for mode != Teuchos::NO_TRANS. More...
 
values_dualv_type::t_host::const_type getValuesViewHost (const RowInfo &rowinfo) const
 Get a const Host view of the locally owned values row myRow, such that rowinfo = getRowInfo(myRow). More...
 
values_dualv_type::t_dev::const_type getValuesViewDevice (const RowInfo &rowinfo) const
 Get a const Device view of the locally owned values row myRow, such that rowinfo = getRowInfo(myRow). More...
 
values_dualv_type::t_host getValuesViewHostNonConst (const RowInfo &rowinfo)
 Get a non-const Host view of the locally owned values row myRow, such that rowinfo = getRowInfo(myRow). More...
 
values_dualv_type::t_dev getValuesViewDeviceNonConst (const RowInfo &rowinfo)
 Get a non-const Device view of the locally owned values row myRow, such that rowinfo = getRowInfo(myRow). More...
 
void swap (CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &matrix)
 Swaps the data from *this with the data and maps from crsMatrix. More...
 
void fillLocalMatrix (const Teuchos::RCP< Teuchos::ParameterList > &params)
 Fill data into the local matrix. More...
 
void fillLocalGraphAndMatrix (const Teuchos::RCP< Teuchos::ParameterList > &params)
 Fill data into the local graph and matrix. More...
 
void checkInternalState () const
 Check that this object's state is sane; throw if it's not. More...
 
virtual size_t constantNumberOfPackets () const
 Whether the implementation's instance promises always to have a constant number of packets per LID (local index), and if so, how many packets per LID there are. More...
 
virtual void doTransfer (const SrcDistObject &src, const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &transfer, const char modeString[], const ReverseOption revOp, const CombineMode CM, const bool restrictedMode)
 Redistribute data across (MPI) processes. More...
 
virtual bool reallocArraysForNumPacketsPerLid (const size_t numExportLIDs, const size_t numImportLIDs)
 Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary. More...
 
void beginTransfer (const SrcDistObject &src, const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &transfer, const char modeString[], const ReverseOption revOp, const CombineMode CM, const bool restrictedMode)
 Implementation detail of doTransfer. More...
 

Static Protected Member Functions

static size_t mergeRowIndicesAndValues (size_t rowLen, local_ordinal_type *cols, impl_scalar_type *vals)
 Merge duplicate row indices in the given row, along with their corresponding values. More...
 

Protected Attributes

Teuchos::RCP< MVimportMV_
 Column Map MultiVector used in apply(). More...
 
Teuchos::RCP< MVexportMV_
 Row Map MultiVector used in apply(). More...
 
Details::EStorageStatus storageStatus_
 Status of the matrix's storage, when not in a fill-complete state. More...
 
bool fillComplete_ = false
 Whether the matrix is fill complete. More...
 
std::map< GlobalOrdinal,
std::pair< Teuchos::Array
< GlobalOrdinal >
, Teuchos::Array< Scalar > > > 
nonlocals_
 Nonlocal data added using insertGlobalValues(). More...
 
(Global) graph pointers

We keep two graph pointers in order to maintain const correctness. myGraph_ is a graph which we create internally. Operations that change the sparsity structure also modify myGraph_. If myGraph_ != null, then staticGraph_ == myGraph_ pointerwise (we set the pointers equal to each other when we create myGraph_). myGraph_ is only null if this CrsMatrix was created using the constructor with a const CrsGraph input argument. In this case, staticGraph_ is set to the input CrsGraph.

Teuchos::RCP< const GraphstaticGraph_
 
Teuchos::RCP< GraphmyGraph_
 

Related Functions

(Note that these are not member functions.)

template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > importAndFillCompleteCrsMatrix (const Teuchos::RCP< const CrsMatrixType > &sourceMatrix, const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &importer, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Nonmember CrsMatrix constructor that fuses Import and fillComplete(). More...
 
template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > importAndFillCompleteCrsMatrix (const Teuchos::RCP< const CrsMatrixType > &sourceMatrix, const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &rowImporter, const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &domainImporter, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
 Nonmember CrsMatrix constructor that fuses Import and fillComplete(). More...
 
template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > exportAndFillCompleteCrsMatrix (const Teuchos::RCP< const CrsMatrixType > &sourceMatrix, const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &exporter, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Nonmember CrsMatrix constructor that fuses Export and fillComplete(). More...
 
template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > exportAndFillCompleteCrsMatrix (const Teuchos::RCP< const CrsMatrixType > &sourceMatrix, const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &rowExporter, const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &domainExporter, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
 Nonmember CrsMatrix constructor that fuses Export and fillComplete(). More...
 
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< CrsMatrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
createCrsMatrix (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const size_t maxNumEntriesPerRow=0, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Create an empty CrsMatrix given a row map and a single integer upper bound on the number of stored entries per row. More...
 
template<class CrsMatrixType >
void removeCrsMatrixZeros (CrsMatrixType &matrix, typename Teuchos::ScalarTraits< typename CrsMatrixType::scalar_type >::magnitudeType const &threshold=Teuchos::ScalarTraits< typename CrsMatrixType::scalar_type >::magnitude(Teuchos::ScalarTraits< typename CrsMatrixType::scalar_type >::zero()))
 Remove zero entries from a matrix. More...
 

Constructors and destructor

template<class S2 , class LO2 , class GO2 , class N2 >
class CrsMatrix
 
template<class S2 , class LO2 , class GO2 , class N2 >
void Details::residual (const Operator< S2, LO2, GO2, N2 > &A, const MultiVector< S2, LO2, GO2, N2 > &X, const MultiVector< S2, LO2, GO2, N2 > &B, MultiVector< S2, LO2, GO2, N2 > &R)
 
template<class MatrixArray , class MultiVectorArray >
void batchedApply (const MatrixArray &Matrices, const typename std::remove_pointer< typename MultiVectorArray::value_type >::type &X, MultiVectorArray &Y, typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type alpha, typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type beta, Teuchos::RCP< Teuchos::ParameterList > params)
 Does multiply matrix apply() calls with a single X vector. More...
 
 CrsMatrix (const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &)=default
 Copy constructor. More...
 
 CrsMatrix (CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &&)=default
 Move constructor. More...
 
CrsMatrixoperator= (const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &)=default
 Copy assignment. More...
 
CrsMatrixoperator= (CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &&)=default
 Move assignment. More...
 
 CrsMatrix (const Teuchos::RCP< const map_type > &rowMap, const size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying the maximum number of entries that any row on the process can take. More...
 
 CrsMatrix (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::ArrayView< const size_t > &numEntPerRowToAlloc, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying (possibly different) number of entries in each row. More...
 
 CrsMatrix (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const size_t maxNumEntPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and fixed number of entries for each row. More...
 
 CrsMatrix (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::ArrayView< const size_t > &numEntPerRowToAlloc, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and number of entries in each row. More...
 
 CrsMatrix (CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &matrix, const Teuchos::RCP< const crs_graph_type > &graph, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying a matrix and a previously constructed graph, presumably a subset of the matrix's graph. This matrix will alias the first N values of the passed-in matrix, where N is the number of entries in the graph. More...
 
 CrsMatrix (const Teuchos::RCP< const crs_graph_type > &graph, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying a previously constructed graph. More...
 
 CrsMatrix (const Teuchos::RCP< const crs_graph_type > &graph, const typename local_matrix_device_type::values_type &values, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying a previously constructed graph and entries array. More...
 
 CrsMatrix (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const typename local_graph_device_type::row_map_type &rowPointers, const typename local_graph_device_type::entries_type::non_const_type &columnIndices, const typename local_matrix_device_type::values_type &values, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and arrays containing the matrix in local indices. In almost all cases the indices must be sorted on input, but if they aren't sorted, "sorted" must be set to false in params. More...
 
 CrsMatrix (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::ArrayRCP< size_t > &rowPointers, const Teuchos::ArrayRCP< LocalOrdinal > &columnIndices, const Teuchos::ArrayRCP< Scalar > &values, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and arrays containing the local matrix. In almost all cases the local matrix must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params. More...
 
 CrsMatrix (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const local_matrix_device_type &lclMatrix, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and a local matrix, which the resulting CrsMatrix views. More...
 
 CrsMatrix (const local_matrix_device_type &lclMatrix, const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap=Teuchos::null, const Teuchos::RCP< const map_type > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column, domain and range Maps, and a local matrix, which the resulting CrsMatrix views. In almost all cases the local matrix must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params. More...
 
 CrsMatrix (const local_matrix_device_type &lclMatrix, const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< const import_type > &importer, const Teuchos::RCP< const export_type > &exporter, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 
 CrsMatrix (const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Teuchos::DataAccess copyOrView)
 Copy constructor, with option to do deep or shallow copy. More...
 
virtual ~CrsMatrix ()=default
 Destructor (virtual for memory safety of derived classes). More...
 

Methods for inserting, modifying, or removing entries

void insertGlobalValues (const GlobalOrdinal globalRow, const Teuchos::ArrayView< const GlobalOrdinal > &cols, const Teuchos::ArrayView< const Scalar > &vals)
 Insert one or more entries into the matrix, using global column indices. More...
 
void insertGlobalValues (const GlobalOrdinal globalRow, const LocalOrdinal numEnt, const Scalar vals[], const GlobalOrdinal inds[])
 Epetra compatibility version of insertGlobalValues (see above) that takes arguments as raw pointers, rather than Teuchos::ArrayView. More...
 
void insertLocalValues (const LocalOrdinal localRow, const Teuchos::ArrayView< const LocalOrdinal > &cols, const Teuchos::ArrayView< const Scalar > &vals, const CombineMode CM=ADD)
 Insert one or more entries into the matrix, using local column indices. More...
 
void insertLocalValues (const LocalOrdinal localRow, const LocalOrdinal numEnt, const Scalar vals[], const LocalOrdinal cols[], const CombineMode CM=ADD)
 Epetra compatibility version of insertLocalValues (see above) that takes arguments as raw pointers, rather than Teuchos::ArrayView. More...
 
local_ordinal_type replaceGlobalValues (const global_ordinal_type globalRow, const Kokkos::View< const global_ordinal_type *, Kokkos::AnonymousSpace > &inputInds, const Kokkos::View< const impl_scalar_type *, Kokkos::AnonymousSpace > &inputVals)
 Replace one or more entries' values, using global indices. More...
 
LocalOrdinal replaceGlobalValues (const GlobalOrdinal globalRow, const Teuchos::ArrayView< const GlobalOrdinal > &cols, const Teuchos::ArrayView< const Scalar > &vals)
 Overload of replaceGlobalValues (see above), that takes Teuchos::ArrayView (host pointers) instead of Kokkos::View. More...
 
LocalOrdinal replaceGlobalValues (const GlobalOrdinal globalRow, const LocalOrdinal numEnt, const Scalar vals[], const GlobalOrdinal cols[])
 Overload of replaceGlobalValues (see above), that takes raw pointers instead of Kokkos::View. More...
 
local_ordinal_type replaceLocalValues (const local_ordinal_type localRow, const Kokkos::View< const local_ordinal_type *, Kokkos::AnonymousSpace > &inputInds, const Kokkos::View< const impl_scalar_type *, Kokkos::AnonymousSpace > &inputVals)
 Replace one or more entries' values, using local row and column indices. More...
 
LocalOrdinal replaceLocalValues (const LocalOrdinal localRow, const Teuchos::ArrayView< const LocalOrdinal > &cols, const Teuchos::ArrayView< const Scalar > &vals)
 Backwards compatibility version of replaceLocalValues (see above), that takes Teuchos::ArrayView (host pointers) instead of Kokkos::View. More...
 
LocalOrdinal replaceLocalValues (const LocalOrdinal localRow, const LocalOrdinal numEnt, const Scalar inputVals[], const LocalOrdinal inputCols[])
 Epetra compatibility version of replaceLocalValues, that takes raw pointers instead of Kokkos::View. More...
 
LocalOrdinal sumIntoGlobalValues (const GlobalOrdinal globalRow, const Teuchos::ArrayView< const GlobalOrdinal > &cols, const Teuchos::ArrayView< const Scalar > &vals, const bool atomic=useAtomicUpdatesByDefault)
 Sum into one or more sparse matrix entries, using global indices. More...
 
LocalOrdinal sumIntoGlobalValues (const GlobalOrdinal globalRow, const LocalOrdinal numEnt, const Scalar vals[], const GlobalOrdinal cols[], const bool atomic=useAtomicUpdatesByDefault)
 Epetra compatibility version of sumIntoGlobalValues (see above), that takes input as raw pointers instead of Kokkos::View. More...
 
local_ordinal_type sumIntoLocalValues (const local_ordinal_type localRow, const Kokkos::View< const local_ordinal_type *, Kokkos::AnonymousSpace > &inputInds, const Kokkos::View< const impl_scalar_type *, Kokkos::AnonymousSpace > &inputVals, const bool atomic=useAtomicUpdatesByDefault)
 Sum into one or more sparse matrix entries, using local row and column indices. More...
 
LocalOrdinal sumIntoLocalValues (const LocalOrdinal localRow, const Teuchos::ArrayView< const LocalOrdinal > &cols, const Teuchos::ArrayView< const Scalar > &vals, const bool atomic=useAtomicUpdatesByDefault)
 Sum into one or more sparse matrix entries, using local row and column indices. More...
 
LocalOrdinal sumIntoLocalValues (const LocalOrdinal localRow, const LocalOrdinal numEnt, const Scalar vals[], const LocalOrdinal cols[], const bool atomic=useAtomicUpdatesByDefault)
 Epetra compatibility version of sumIntoLocalValues (see above) that takes raw pointers instead of Kokkos::View. More...
 
template<class LocalIndicesViewType , class ImplScalarViewType , class BinaryFunction >
LocalOrdinal transformLocalValues (const LocalOrdinal lclRow, const typename UnmanagedView< LocalIndicesViewType >::type &inputInds, const typename UnmanagedView< ImplScalarViewType >::type &inputVals, BinaryFunction f, const bool atomic=useAtomicUpdatesByDefault)
 Transform CrsMatrix entries in place, using local indices to select the entries in the row to transform. More...
 
template<class BinaryFunction , class InputMemorySpace >
LocalOrdinal transformGlobalValues (const GlobalOrdinal gblRow, const Kokkos::View< const GlobalOrdinal *, InputMemorySpace, Kokkos::MemoryUnmanaged > &inputInds, const Kokkos::View< const impl_scalar_type *, InputMemorySpace, Kokkos::MemoryUnmanaged > &inputVals, BinaryFunction f, const bool atomic=useAtomicUpdatesByDefault)
 Transform CrsMatrix entries in place, using global indices to select the entries in the row to transform. More...
 
void setAllToScalar (const Scalar &alpha)
 Set all matrix entries equal to alpha. More...
 
void scale (const Scalar &alpha)
 Scale the matrix's values: this := alpha*this. More...
 
void setAllValues (const typename local_graph_device_type::row_map_type &ptr, const typename local_graph_device_type::entries_type::non_const_type &ind, const typename local_matrix_device_type::values_type &val)
 Set the local matrix using three (compressed sparse row) arrays. More...
 
void setAllValues (const local_matrix_device_type &localMatrix)
 Set the local matrix using an existing local matrix. More...
 
void setAllValues (const Teuchos::ArrayRCP< size_t > &ptr, const Teuchos::ArrayRCP< LocalOrdinal > &ind, const Teuchos::ArrayRCP< Scalar > &val)
 Set the local matrix using three (compressed sparse row) arrays. More...
 
row_ptrs_host_view_type getLocalRowPtrsHost () const
 Get a host view of the CRS packed row pointers. More...
 
row_ptrs_device_view_type getLocalRowPtrsDevice () const
 Get a device view of the CRS packed row pointers. More...
 
local_inds_host_view_type getLocalIndicesHost () const
 Get a host view of the CRS packed column indicies. More...
 
local_inds_device_view_type getLocalIndicesDevice () const
 Get a device_view of the CRS packed column indicies. More...
 
virtual LocalOrdinal replaceGlobalValuesImpl (impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const GlobalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts)
 Implementation detail of replaceGlobalValues. More...
 
virtual LocalOrdinal replaceLocalValuesImpl (impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const LocalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts)
 Implementation detail of replaceLocalValues. More...
 
virtual LocalOrdinal sumIntoGlobalValuesImpl (impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const GlobalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts, const bool atomic=useAtomicUpdatesByDefault)
 Implementation detail of sumIntoGlobalValues. More...
 
virtual LocalOrdinal sumIntoLocalValuesImpl (impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const LocalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts, const bool atomic=useAtomicUpdatesByDefault)
 Implementation detail of sumIntoLocalValues. More...
 

Methods implementing RowMatrix

using values_dualv_type = Kokkos::DualView< impl_scalar_type *, device_type >
 
using values_wdv_type = Details::WrappedDualView< values_dualv_type >
 
using ordinal_rowptrs_type = typename local_multiply_op_type::ordinal_view_type
 
values_wdv_type valuesUnpacked_wdv
 
values_wdv_type valuesPacked_wdv
 
ordinal_rowptrs_type ordinalRowptrs
 local_ordinal typed version of local matrix's rowptrs. This allows the LocalCrsMatrixOperator to have rowptrs and entries be the same type, so cuSPARSE SpMV (including merge-path) can be used for apply. This is allocated and populated lazily in getLocalMultiplyOperator(), only if all 4 conditions are met: More...
 
Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm () const override
 The communicator over which the matrix is distributed. More...
 
Teuchos::RCP< const map_typegetRowMap () const override
 The Map that describes the row distribution in this matrix. More...
 
Teuchos::RCP< const map_typegetColMap () const override
 The Map that describes the column distribution in this matrix. More...
 
Teuchos::RCP< const RowGraph
< LocalOrdinal, GlobalOrdinal,
Node > > 
getGraph () const override
 This matrix's graph, as a RowGraph. More...
 
Teuchos::RCP< const
crs_graph_type
getCrsGraph () const
 This matrix's graph, as a CrsGraph. More...
 
local_matrix_device_type getLocalMatrixDevice () const
 The local sparse matrix. More...
 
local_matrix_host_type getLocalMatrixHost () const
 
std::shared_ptr
< local_multiply_op_type
getLocalMultiplyOperator () const
 The local sparse matrix operator (a wrapper of getLocalMatrixDevice() that supports local matrix-vector multiply) More...
 
global_size_t getGlobalNumRows () const override
 Number of global elements in the row map of this matrix. More...
 
global_size_t getGlobalNumCols () const override
 The number of global columns in the matrix. More...
 
size_t getLocalNumRows () const override
 The number of matrix rows owned by the calling process. More...
 
size_t getLocalNumCols () const override
 The number of columns connected to the locally owned rows of this matrix. More...
 
GlobalOrdinal getIndexBase () const override
 The index base for global indices for this matrix. More...
 
global_size_t getGlobalNumEntries () const override
 The global number of entries in this matrix. More...
 
size_t getLocalNumEntries () const override
 The local number of entries in this matrix. More...
 
size_t getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const override
 Number of entries in the sparse matrix in the given global row, on the calling (MPI) process. More...
 
size_t getNumEntriesInLocalRow (local_ordinal_type localRow) const override
 Number of entries in the sparse matrix in the given local row, on the calling (MPI) process. More...
 
size_t getGlobalMaxNumRowEntries () const override
 Maximum number of entries in any row of the matrix, over all processes in the matrix's communicator. More...
 
size_t getLocalMaxNumRowEntries () const override
 Maximum number of entries in any row of the matrix, on this process. More...
 
virtual LocalOrdinal getBlockSize () const override
 The number of degrees of freedom per mesh point. More...
 
bool hasColMap () const override
 Whether the matrix has a well-defined column Map. More...
 
bool isLocallyIndexed () const override
 Whether the matrix is locally indexed on the calling process. More...
 
bool isGloballyIndexed () const override
 Whether the matrix is globally indexed on the calling process. More...
 
bool isFillComplete () const override
 Whether the matrix is fill complete. More...
 
bool isFillActive () const
 Whether the matrix is not fill complete. More...
 
bool isStorageOptimized () const
 Returns true if storage has been optimized. More...
 
bool isStaticGraph () const
 Indicates that the graph is static, so that new entries cannot be added to this matrix. More...
 
mag_type getFrobeniusNorm () const override
 Compute and return the Frobenius norm of the matrix. More...
 
virtual bool supportsRowViews () const override
 Return true if getLocalRowView() and getGlobalRowView() are valid for this object. More...
 
void getGlobalRowCopy (GlobalOrdinal GlobalRow, nonconst_global_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const override
 Fill given arrays with a deep copy of the locally owned entries of the matrix in a given row, using global column indices. More...
 
void getLocalRowCopy (LocalOrdinal LocalRow, nonconst_local_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const override
 Fill given arrays with a deep copy of the locally owned entries of the matrix in a given row, using local column indices. More...
 
void getGlobalRowView (GlobalOrdinal GlobalRow, global_inds_host_view_type &indices, values_host_view_type &values) const override
 Get a constant, nonpersisting view of a row of this matrix, using global row and column indices. More...
 
void getLocalRowView (LocalOrdinal LocalRow, local_inds_host_view_type &indices, values_host_view_type &values) const override
 Get a constant view of a row of this matrix, using local row and column indices. More...
 
void getLocalDiagCopy (Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const override
 Get a constant, nonpersisting view of a row of this matrix, using local row and column indices, with raw pointers. More...
 
void getLocalDiagOffsets (Teuchos::ArrayRCP< size_t > &offsets) const
 Get offsets of the diagonal entries in the matrix. More...
 
void getLocalDiagCopy (Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Kokkos::View< const size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
 Variant of getLocalDiagCopy() that uses precomputed offsets. More...
 
void getLocalDiagCopy (Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView< const size_t > &offsets) const
 Variant of getLocalDiagCopy() that uses precomputed offsets. More...
 
void leftScale (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) override
 Scale the matrix on the left with the given Vector. More...
 
void rightScale (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) override
 Scale the matrix on the right with the given Vector. More...
 

Implementation of DistObject interface

typedef DistObject< Scalar,
LocalOrdinal, GlobalOrdinal,
Node >::buffer_device_type 
buffer_device_type
 Kokkos::Device specialization for communication buffers. More...
 
virtual bool checkSizes (const SrcDistObject &source) override
 Compare the source and target (this) objects for compatibility. More...
 
void applyCrsPadding (const typename crs_graph_type::padding_type &padding, const bool verbose)
 
void unpackAndCombine (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< char *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode CM) override
 Unpack the imported column indices and values, and combine into matrix. More...
 
void packNew (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< char *, buffer_device_type > &exports, const Kokkos::DualView< size_t *, buffer_device_type > &numPacketsPerLID, size_t &constantNumPackets) const
 Pack this object's data for an Import or Export. More...
 
virtual void copyAndPermute (const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM) override
 
virtual void packAndPrepare (const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< char *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets) override
 

Methods implemented by subclasses and used by doTransfer().

The doTransfer() method uses the subclass' implementations of these methods to implement data transfer. Subclasses of DistObject must implement these methods. This is an instance of the Template Method Pattern. ("Template" here doesn't mean "C++ template"; it means "pattern with holes that are filled in by the subclass' method implementations.")

Teuchos::RCP< const map_typemap_
 The Map over which this object is distributed. More...
 
Kokkos::DualView< packet_type
*, buffer_device_type
imports_
 Buffer into which packed data are imported (received from other processes). More...
 
Kokkos::DualView< size_t
*, buffer_device_type
numImportPacketsPerLID_
 Number of packets to receive for each receive operation. More...
 
Kokkos::DualView< packet_type
*, buffer_device_type
exports_
 Buffer from which packed data are exported (sent to other processes). More...
 
Kokkos::DualView< size_t
*, buffer_device_type
numExportPacketsPerLID_
 Number of packets to send for each send operation. More...
 
virtual void copyAndPermute (const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM)
 Perform copies and permutations that are local to the calling (MPI) process. More...
 
virtual void copyAndPermute (const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM, const execution_space &space)
 Same as copyAndPermute, but do operations in space. More...
 
virtual void packAndPrepare (const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets)
 Pack data and metadata for communication (sends). More...
 
virtual void packAndPrepare (const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets, const execution_space &space)
 Same as packAndPrepare, but in an execution space instance. More...
 
virtual void unpackAndCombine (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode)
 Perform any unpacking and combining after communication. More...
 
virtual void unpackAndCombine (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode, const execution_space &space)
 
std::unique_ptr< std::string > createPrefix (const char className[], const char methodName[]) const
 
virtual bool reallocImportsIfNeeded (const size_t newSize, const bool verbose, const std::string *prefix, const bool remoteLIDsContiguous=false, const CombineMode CM=INSERT)
 Reallocate imports_ if needed. More...
 

Detailed Description

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
class Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >

Sparse matrix that presents a row-oriented interface that lets users read or modify entries.

Template Parameters
ScalarThe type of the numerical entries of the matrix. (You can use real-valued or complex-valued types here, unlike in Epetra, where the scalar type is always double.)
LocalOrdinalThe type of local indices. See the documentation of Map for requirements.
GlobalOrdinalThe type of global indices. See the documentation of Map for requirements.
NodeThe Kokkos Node type. See the documentation of Map for requirements.

This class implements a distributed-memory parallel sparse matrix, and provides sparse matrix-vector multiply (including transpose) and sparse triangular solve operations. It provides access by rows to the elements of the matrix, as if the local data were stored in compressed sparse row format. (Implementations are not required to store the data in this way internally.) This class has an interface like that of Epetra_CrsMatrix, but also allows insertion of data into nonowned rows, much like Epetra_FECrsMatrix.

Prerequisites

Before reading the rest of this documentation, it helps to know something about the Teuchos memory management classes, in particular Teuchos::RCP, Teuchos::ArrayRCP, and Teuchos::ArrayView. You should also know a little bit about MPI (the Message Passing Interface for distributed-memory programming). You won't have to use MPI directly to use CrsMatrix, but it helps to be familiar with the general idea of distributed storage of data over a communicator. Finally, you should read the documentation of Map and MultiVector.

Local and global indices

The distinction between local and global indices might confuse new Tpetra users. Please refer to the documentation of Map for a detailed explanation. This is important because many of CrsMatrix's methods for adding, modifying, or accessing entries come in versions that take either local or global indices. The matrix itself may store indices either as local or global, and the same matrix may use global indices or local indices at different points in its life. You should only use the method version corresponding to the current state of the matrix. For example, getGlobalRowView() returns a view to the indices represented as global; it is incorrect to call this method if the matrix is storing indices as local. Call isGloballyIndexed() or isLocallyIndexed() to find out whether the matrix currently stores indices as local or global.

It may also help to read CrsGraph's documentation.

Insertion into nonowned rows

All methods (except for insertGlobalValues() and sumIntoGlobalValues(); see below) that work with global indices only allow operations on indices owned by the calling process. For example, methods that take a global row index expect that row to be owned by the calling process. Access to nonowned rows, that is, rows not owned by the calling process, requires performing an explicit communication via the Import / Export capabilities of the CrsMatrix object. See the documentation of DistObject for more details.

The methods insertGlobalValues() and sumIntoGlobalValues() are exceptions to this rule. They both allows you to add data to nonowned rows. These data are stored locally and communicated to the appropriate process on the next call to globalAssemble() or fillComplete(). This means that CrsMatrix provides the same nonowned insertion functionality that Epetra provides via Epetra_FECrsMatrix.

Note for developers on DistObject

DistObject only takes a single Map as input to its constructor. MultiVector is an example of a subclass for which a single Map suffices to describe its data distribution. In that case, DistObject's getMap() method obviously must return that Map. CrsMatrix is an example of a subclass that requires two Map objects: a row Map and a column Map. For CrsMatrix, getMap() returns the row Map. This means that doTransfer() (which CrsMatrix does not override) uses the row Map objects of the source and target CrsMatrix objects. CrsMatrix in turn uses its column Map (if it has one) to "filter" incoming sparse matrix entries whose column indices are not in that process' column Map. This means that CrsMatrix may perform extra communication, though the Import and Export operations are still correct.

This is necessary if the CrsMatrix does not yet have a column Map. Other processes might have added new entries to the matrix; the calling process has to see them in order to accept them. However, the CrsMatrix may already have a column Map, for example, if it was created with the constructor that takes both a row and a column Map, or if it is fill complete (which creates the column Map if the matrix does not yet have one). In this case, it could be possible to "filter" on the sender (instead of on the receiver, as CrsMatrix currently does) and avoid sending data corresponding to columns that the receiver does not own. Doing this would require revising the Import or Export object (instead of the incoming data) using the column Map, to remove global indices and their target process ranks from the send lists if the target process does not own those columns, and to remove global indices and their source process ranks from the receive lists if the calling process does not own those columns. (Abstractly, this is a kind of set difference between an Import or Export object for the row Maps, and the Import resp. Export object for the column Maps.) This could be done separate from DistObject, by creating a new "filtered" Import or Export object, that keeps the same source and target Map objects but has a different communication plan. We have not yet implemented this optimization.

Definition at line 426 of file Tpetra_CrsMatrix_decl.hpp.

Member Typedef Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scalar_type = Scalar

The type of each entry in the matrix.

Definition at line 444 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type = LocalOrdinal

The type of each local index in the matrix.

Definition at line 446 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type = GlobalOrdinal

The type of each global index in the matrix.

Definition at line 448 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::device_type = typename Node::device_type

The Kokkos device type.

Definition at line 450 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::execution_space = typename device_type::execution_space

The Kokkos execution space.

Definition at line 452 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::memory_space = typename device_type::memory_space

The Kokkos memory space.

Definition at line 454 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::node_type = Node

This class' Kokkos Node type.

This is a leftover that will be deprecated and removed. See e.g., GitHub Issue #57.

Definition at line 460 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::map_type = Map<LocalOrdinal, GlobalOrdinal, Node>

The Map specialization suitable for this CrsMatrix specialization.

Definition at line 463 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::import_type = Import<LocalOrdinal, GlobalOrdinal, Node>

The Import specialization suitable for this CrsMatrix specialization.

Definition at line 466 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::export_type = Export<LocalOrdinal, GlobalOrdinal, Node>

The Export specialization suitable for this CrsMatrix specialization.

Definition at line 469 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::row_matrix_type = RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>

The RowMatrix representing the base class of CrsMatrix.

Definition at line 472 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::impl_scalar_type = typename row_matrix_type::impl_scalar_type

The type used internally in place of Scalar.

Some Scalar types might not work with Kokkos on all execution spaces, due to missing CUDA device macros or volatile overloads. The C++ standard type std::complex<T> has this problem. To fix this, we replace std::complex<T> values internally with the (usually) bitwise identical type Kokkos::complex<T>. The latter is the impl_scalar_type corresponding to Scalar = std::complex.

Definition at line 483 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type = typename Kokkos::ArithTraits<impl_scalar_type>::mag_type

Type of a norm result.

This is usually the same as the type of the magnitude (absolute value) of Scalar, but may differ for certain Scalar types.

Definition at line 489 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::crs_graph_type = CrsGraph<LocalOrdinal, GlobalOrdinal, Node>

The CrsGraph specialization suitable for this CrsMatrix specialization.

Definition at line 492 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_graph_device_type = typename crs_graph_type::local_graph_device_type

The part of the sparse matrix's graph on each MPI process.

Definition at line 495 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_device_type = KokkosSparse::CrsMatrix<impl_scalar_type, local_ordinal_type, device_type, void, typename local_graph_device_type::size_type>

The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI process.

Definition at line 505 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_multiply_op_type = LocalCrsMatrixOperator<scalar_type, scalar_type, device_type>

The type of the local matrix-vector operator (a wrapper of KokkosSparse::CrsMatrix )

Definition at line 514 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::buffer_device_type

Kokkos::Device specialization for communication buffers.

See #1088 for why this is not just device_type::device_type.

Definition at line 2962 of file Tpetra_CrsMatrix_decl.hpp.

using Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::packet_type = typename ::Kokkos::ArithTraits<char >::val_type
inherited

The type of each datum being sent or received in an Import or Export.

Note that this type does not always correspond to the Scalar template parameter of subclasses.

Definition at line 334 of file Tpetra_DistObject_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  )
default

Copy constructor.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &&  )
default

Move constructor.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const Teuchos::RCP< const map_type > &  rowMap,
const size_t  maxNumEntriesPerRow,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying the maximum number of entries that any row on the process can take.

Parameters
rowMap[in] Distribution of rows of the matrix.
maxNumEntriesPerRow[in] Maximum number of matrix entries per row. This is a strict upper bound. It may differ on different processes.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 143 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::ArrayView< const size_t > &  numEntPerRowToAlloc,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying (possibly different) number of entries in each row.

Parameters
rowMap[in] Distribution of rows of the matrix.
numEntPerRowToAlloc[in] Maximum number of matrix entries to allocate for each row. This is a strict upper bound. It may differ on different processes.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 172 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const size_t  maxNumEntPerRow,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and fixed number of entries for each row.

The column Map will be used to filter any matrix entries inserted using insertLocalValues() or insertGlobalValues().

Parameters
rowMap[in] Distribution of rows of the matrix.
colMap[in] Distribution of columns of the matrix. See replaceColMap() for the requirements.
maxNumEntPerRow[in] Maximum number of matrix entries per row. This is a strict upper bound. It may differ on different processes.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 203 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const Teuchos::ArrayView< const size_t > &  numEntPerRowToAlloc,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and number of entries in each row.

The column Map will be used to filter any matrix indices inserted using insertLocalValues() or insertGlobalValues().

Parameters
rowMap[in] Distribution of rows of the matrix.
colMap[in] Distribution of columns of the matrix. See replaceColMap() for the requirements.
numEntPerRowToAlloc[in] Maximum number of matrix entries to allocate for each row. This is a strict upper bound. It may differ on different processes.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 246 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  matrix,
const Teuchos::RCP< const crs_graph_type > &  graph,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)
explicit

Constructor specifying a matrix and a previously constructed graph, presumably a subset of the matrix's graph. This matrix will alias the first N values of the passed-in matrix, where N is the number of entries in the graph.

Calling this constructor fixes the graph structure of the sparse matrix. We say in this case that the matrix has a "static graph." If you create a CrsMatrix with this constructor, you are not allowed to insert new entries into the matrix, but you are allowed to change values in the matrix.

The given graph must be fill complete. Note that calling resumeFill() on the graph makes it not fill complete, even if you had previously called fillComplete() on the graph. In that case, you must call fillComplete() on the graph again before invoking this CrsMatrix constructor.

This constructor is marked explicit so that you can't create a CrsMatrix by accident when passing a CrsGraph into a function that takes a CrsMatrix.

Parameters
matrix[in] The existing matrix whose values this one will alias.
graph[in] The graph structure of the sparse matrix. The graph must be fill complete.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 337 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const Teuchos::RCP< const crs_graph_type > &  graph,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)
explicit

Constructor specifying a previously constructed graph.

Calling this constructor fixes the graph structure of the sparse matrix. We say in this case that the matrix has a "static graph." If you create a CrsMatrix with this constructor, you are not allowed to insert new entries into the matrix, but you are allowed to change values in the matrix.

The given graph must be fill complete. Note that calling resumeFill() on the graph makes it not fill complete, even if you had previously called fillComplete() on the graph. In that case, you must call fillComplete() on the graph again before invoking this CrsMatrix constructor.

This constructor is marked explicit so that you can't create a CrsMatrix by accident when passing a CrsGraph into a function that takes a CrsMatrix.

Parameters
graph[in] The graph structure of the sparse matrix. The graph must be fill complete.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 279 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const Teuchos::RCP< const crs_graph_type > &  graph,
const typename local_matrix_device_type::values_type &  values,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)
explicit

Constructor specifying a previously constructed graph and entries array.

Calling this constructor fixes the graph structure of the sparse matrix. We say in this case that the matrix has a "static graph." If you create a CrsMatrix with this constructor, you are not allowed to insert new entries into the matrix, but you are allowed to change values in the matrix.

The given graph must be fill complete. Note that calling resumeFill() on the graph makes it not fill complete, even if you had previously called fillComplete() on the graph. In that case, you must call fillComplete() on the graph again before invoking this CrsMatrix constructor.

This constructor is marked explicit so that you can't create a CrsMatrix by accident when passing a CrsGraph into a function that takes a CrsMatrix.

Parameters
graph[in] The graph structure of the sparse matrix. The graph must be fill complete.
values[in] The local entries in the matrix, as in a CSR "vals" array. The length of this vector should be equal to the number of unknowns in the matrix.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 369 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const typename local_graph_device_type::row_map_type &  rowPointers,
const typename local_graph_device_type::entries_type::non_const_type &  columnIndices,
const typename local_matrix_device_type::values_type &  values,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and arrays containing the matrix in local indices. In almost all cases the indices must be sorted on input, but if they aren't sorted, "sorted" must be set to false in params.

Parameters
rowMap[in] Distribution of rows of the matrix.
colMap[in] Distribution of columns of the matrix. See replaceColMap() for the requirements.
rowPointers[in] The beginning of each row in the matrix, as in a CSR "rowptr" array. The length of this vector should be equal to the number of rows in the graph, plus one. This last entry should store the nunber of nonzeros in the matrix.
columnIndices[in] The local indices of the columns, as in a CSR "colind" array. The length of this vector should be equal to the number of unknowns in the matrix.
values[in] The local entries in the matrix, as in a CSR "vals" array. The length of this vector should be equal to the number of unknowns in the matrix.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 409 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const Teuchos::ArrayRCP< size_t > &  rowPointers,
const Teuchos::ArrayRCP< LocalOrdinal > &  columnIndices,
const Teuchos::ArrayRCP< Scalar > &  values,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and arrays containing the local matrix. In almost all cases the local matrix must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params.

Parameters
rowMap[in] Distribution of rows of the matrix.
colMap[in] Distribution of columns of the matrix. See replaceColMap() for the requirements.
rowPointers[in] The beginning of each row in the matrix, as in a CSR "rowptr" array. The length of this vector should be equal to the number of rows in the graph, plus one. This last entry should store the nunber of nonzeros in the matrix.
columnIndices[in] The local indices of the columns, as in a CSR "colind" array. The length of this vector should be equal to the number of unknowns in the matrix.
values[in] The local entries in the matrix, as in a CSR "vals" array. The length of this vector should be equal to the number of unknowns in the matrix.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 520 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const local_matrix_device_type lclMatrix,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and a local matrix, which the resulting CrsMatrix views.

Unlike most other CrsMatrix constructors, successful completion of this constructor will result in a fill-complete matrix.

Parameters
rowMap[in] Distribution of rows of the matrix.
colMap[in] Distribution of columns of the matrix. See replaceColMap() for the requirements.
lclMatrix[in] A local CrsMatrix containing all local matrix values as well as a local graph. The graph's local row indices must come from the specified row Map, and its local column indices must come from the specified column Map.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 588 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const local_matrix_device_type lclMatrix,
const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const Teuchos::RCP< const map_type > &  domainMap = Teuchos::null,
const Teuchos::RCP< const map_type > &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column, domain and range Maps, and a local matrix, which the resulting CrsMatrix views. In almost all cases the local matrix must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params.

Unlike most other CrsMatrix constructors, successful completion of this constructor will result in a fill-complete matrix.

Parameters
rowMap[in] Distribution of rows of the matrix.
colMap[in] Distribution of columns of the matrix. See replaceColMap() for the requirements.
domainMap[in] The matrix's domain Map. MUST be one to one!
rangeMap[in] The matrix's range Map. MUST be one to one! May be, but need not be, the same as the domain Map.
lclMatrix[in] A local CrsMatrix containing all local matrix values as well as a local graph. The graph's local row indices must come from the specified row Map, and its local column indices must come from the specified column Map.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 639 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const local_matrix_device_type lclMatrix,
const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const Teuchos::RCP< const import_type > &  importer,
const Teuchos::RCP< const export_type > &  exporter,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Create a fill-complete CrsMatrix from all the things it needs.

Parameters
lclMatrix[in] In almost all cases the local matrix must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params.

Definition at line 693 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CrsMatrix ( const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Teuchos::DataAccess  copyOrView 
)

Copy constructor, with option to do deep or shallow copy.

Definition at line 750 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::~CrsMatrix ( )
virtualdefault

Destructor (virtual for memory safety of derived classes).

Note
To Tpetra developers: See the C++ Core Guidelines C.21 ("If you define or <tt>=delete</tt> any default operation, define or <tt>=delete</tt> them all"), in particular the AbstractBase example, for why this destructor declaration implies that we need the above four =default declarations for copy construction, move construction, copy assignment, and move assignment.

Member Function Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
CrsMatrix& Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::operator= ( const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  )
default

Copy assignment.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
CrsMatrix& Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::operator= ( CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &&  )
default

Move assignment.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal, class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::insertGlobalValues ( const GlobalOrdinal  globalRow,
const Teuchos::ArrayView< const GlobalOrdinal > &  cols,
const Teuchos::ArrayView< const Scalar > &  vals 
)

Insert one or more entries into the matrix, using global column indices.

Parameters
globalRow[in] Global index of the row into which to insert the entries.
cols[in] Global indices of the columns into which to insert the entries.
vals[in] Values to insert into the above columns.

For all k in 0, ..., col.size()-1, insert the value values[k] into entry (globalRow, cols[k]) of the matrix. If that entry already exists, add the new value to the old value.

This is a local operation. It does not communicate (using MPI). If row globalRow is owned by the calling process, the entries will be inserted immediately. Otherwise, if that row is not owned by the calling process, then the entries will be stored locally for now, and only communicated to the process that owns the row when either fillComplete() or globalAssemble() is called. If that process already has an entry, the incoming value will be added to the old value, just as if it were inserted on the owning process. If the matrix has a column Map (hasColMap() == true), and if globalRow is owned by process p, then it is forbidden to insert column indices that are not in the column Map on process p. Tpetra will test the input column indices to ensure this is the case, but if globalRow is not owned by the calling process, the test will be deferred until the next call to globalAssemble() or fillComplete().

Warning
The behavior described in the above paragraph differs from that of Epetra. If the matrix has a column Map, Epetra_CrsMatrix "filters" column indices not in the column Map. Many users found this confusing, so we changed it so that nonowned column indices are forbidden.

It is legal to call this method whether the matrix's column indices are globally or locally indexed. If the matrix's column indices are locally indexed (isLocallyIndexed() == true), then this method will convert the input global column indices to local column indices.

For better performance when filling entries into a sparse matrix, consider the following tips:

  1. Use local indices (e.g., insertLocalValues()) if you know the column Map in advance. Converting global indices to local indices is expensive. Of course, if you don't know the column Map in advance, you must use global indices.
  2. When invoking the CrsMatrix constructor, give the best possible upper bounds on the number of entries in each row of the matrix. This will avoid expensive reallocation if your bound was not large enough.
  3. If you plan to reuse a matrix's graph structure, but change its values, in repeated fillComplete() / resumeFill() cycles, you can get the best performance by creating the matrix with a const CrsGraph. Do this by using the CrsMatrix constructor that accepts an RCP of a const CrsGraph. If you do this, you must use the "replace" or "sumInto" methods to change the values of the matrix; you may not use insertGlobalValues() or insertLocalValues().

Definition at line 2007 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::insertGlobalValues ( const GlobalOrdinal  globalRow,
const LocalOrdinal  numEnt,
const Scalar  vals[],
const GlobalOrdinal  inds[] 
)

Epetra compatibility version of insertGlobalValues (see above) that takes arguments as raw pointers, rather than Teuchos::ArrayView.

Arguments are the same and in the same order as Epetra_CrsMatrix::InsertGlobalValues.

Parameters
globalRow[in] Global index of the row into which to insert the entries.
numEnt[in] Number of entries to insert; number of valid entries in vals and inds.
vals[in] Values to insert.
inds[in] Global indices of the columns into which to insert the entries.

Definition at line 2119 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::insertLocalValues ( const LocalOrdinal  localRow,
const Teuchos::ArrayView< const LocalOrdinal > &  cols,
const Teuchos::ArrayView< const Scalar > &  vals,
const CombineMode  CM = ADD 
)

Insert one or more entries into the matrix, using local column indices.

Parameters
localRow[in] Local index of the row into which to insert the entries. It must be owned by the row Map on the calling process.
cols[in] Local indices of the columns into which to insert the entries. All of the column indices must be owned by the column Map on the calling process.
vals[in] Values to insert into the above columns.
CM[in] How values should be inserted. Valid options are: ADD (default) inserts values that are not yet in the matrix graph, and sums values that are already present. INSERT inserts values that are not yet in the matrix graph, and replaces values that are already present.

For all k in 0, ..., cols.size()-1, insert the value values[k] into entry (globalRow, cols[k]) of the matrix. If that entry already exists, add the new value to the old value, if CM=ADD, otherwise replace the old value.

In order to call this method, the matrix must be locally indexed, and it must have a column Map.

For better performance when filling entries into a sparse matrix, consider the following tips:

  1. When invoking the CrsMatrix constructor, give the best possible upper bounds on the number of entries in each row of the matrix. This will avoid expensive reallocation if your bound was not large enough.
  2. If you plan to reuse a matrix's graph structure, but change its values, in repeated fillComplete() / resumeFill() cycles, you can get the best performance by creating the matrix with a const CrsGraph. Do this by using the CrsMatrix constructor that accepts an RCP of a const CrsGraph. If you do this, you must use the "replace" or "sumInto" methods to change the values of the matrix; you may not use insertGlobalValues() or insertLocalValues().

Definition at line 1775 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::insertLocalValues ( const LocalOrdinal  localRow,
const LocalOrdinal  numEnt,
const Scalar  vals[],
const LocalOrdinal  cols[],
const CombineMode  CM = ADD 
)

Epetra compatibility version of insertLocalValues (see above) that takes arguments as raw pointers, rather than Teuchos::ArrayView.

Arguments are the same and in the same order as Epetra_CrsMatrix::InsertMyValues.

Parameters
localRow[in] Local index of the row into which to insert the entries.
numEnt[in] Number of entries to insert; number of valid entries in vals and cols.
vals[in] Values to insert.
cols[in] Global indices of the columns into which to insert the entries.
CM[in] How values should be inserted. Valid options are: ADD (default) inserts values that are not yet in the matrix graph, and sums values that are already present. INSERT inserts values that are not yet in the matrix graph, and replaces values that are already present.

Definition at line 1872 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal, class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceGlobalValuesImpl ( impl_scalar_type  rowVals[],
const crs_graph_type graph,
const RowInfo rowInfo,
const GlobalOrdinal  inds[],
const impl_scalar_type  newVals[],
const LocalOrdinal  numElts 
)
protectedvirtual

Implementation detail of replaceGlobalValues.

Parameters
rowVals[in/out] On input: Values of the row of the sparse matrix to modify. On output: The modified values.
graph[in] The matrix's graph.
rowInfo[in] Result of graph.getRowInfo on the index of the local row of the matrix to modify.
inds[in] Global column indices of that row to modify.
newVals[in] For each k, replace the value in rowVals corresponding to local column index inds[k] with newVals[k].

Definition at line 2449 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceGlobalValues ( const global_ordinal_type  globalRow,
const Kokkos::View< const global_ordinal_type *, Kokkos::AnonymousSpace > &  inputInds,
const Kokkos::View< const impl_scalar_type *, Kokkos::AnonymousSpace > &  inputVals 
)

Replace one or more entries' values, using global indices.

Parameters
globalRow[in] Global index of the row in which to replace the entries. This row must be owned by the calling process.
inputInds[in] Kokkos::View of the global indices of the columns in which to replace the entries.
inputVals[in] Kokkos::View of the values to use for replacing the entries.

For all k in 0, ..., inputInds.extent(0)-1, replace the value at entry (globalRow, inputInds(k)) of the matrix with inputVals(k). That entry must exist in the matrix already.

If (globalRow, inputInds(k)) corresponds to an entry that is duplicated in this matrix row (likely because it was inserted more than once and fillComplete() has not been called in the interim), the behavior of this method is not defined.

Returns
The number of indices for which values were actually replaced; the number of "correct" indices.

If the returned value N satisfies

0 <= N < inputInds.extent(0),

then inputInds.extent(0) - N of the entries of cols are not valid global column indices. If the returned value is Teuchos::OrdinalTraits<LocalOrdinal>::invalid(), then at least one of the following is true:

  • ! isFillActive ()
  • inputInds.extent (0) != inputVals.extent (0)

Definition at line 2517 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal, class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceGlobalValues ( const GlobalOrdinal  globalRow,
const Teuchos::ArrayView< const GlobalOrdinal > &  cols,
const Teuchos::ArrayView< const Scalar > &  vals 
)

Overload of replaceGlobalValues (see above), that takes Teuchos::ArrayView (host pointers) instead of Kokkos::View.

Definition at line 2468 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceGlobalValues ( const GlobalOrdinal  globalRow,
const LocalOrdinal  numEnt,
const Scalar  vals[],
const GlobalOrdinal  cols[] 
)

Overload of replaceGlobalValues (see above), that takes raw pointers instead of Kokkos::View.

This version of the method takes the same arguments in the same order as Epetra_CrsMatrix::ReplaceGlobalValues.

Parameters
globalRow[in] Global index of the row in which to replace the entries. This row must be owned by the calling process.
numEnt[in] Number of entries to replace; number of valid entries in vals and cols.
vals[in] Values to use for replacing the entries.
cols[in] Global indices of the columns in which to replace the entries.

Definition at line 2486 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal , class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceLocalValuesImpl ( impl_scalar_type  rowVals[],
const crs_graph_type graph,
const RowInfo rowInfo,
const LocalOrdinal  inds[],
const impl_scalar_type  newVals[],
const LocalOrdinal  numElts 
)
protectedvirtual

Implementation detail of replaceLocalValues.

Parameters
rowVals[in/out] On input: Values of the row of the sparse matrix to modify. On output: The modified values.
graph[in] The matrix's graph.
rowInfo[in] Result of graph.getRowInfo on the index of the local row of the matrix to modify.
inds[in] Local column indices of that row to modify.
newVals[in] For each k, replace the value in rowVals corresponding to local column index inds[k] with newVals[k].

Definition at line 2311 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceLocalValues ( const local_ordinal_type  localRow,
const Kokkos::View< const local_ordinal_type *, Kokkos::AnonymousSpace > &  inputInds,
const Kokkos::View< const impl_scalar_type *, Kokkos::AnonymousSpace > &  inputVals 
)

Replace one or more entries' values, using local row and column indices.

Parameters
localRow[in] local index of the row in which to replace the entries. This row must be owned by the calling process.
cols[in] Local indices of the columns in which to replace the entries.
vals[in] Values to use for replacing the entries.

For local row index localRow and local column indices cols, do A(localRow, cols(k)) = vals(k). The row index and column indices must be valid on the calling process, and all matrix entries A(localRow, cols(k)) must already exist. (This method does not change the matrix's structure.) If the row index is valid, any invalid column indices are ignored, but counted in the return value.

Returns
The number of indices for which values were actually replaced; the number of "correct" indices.

If the returned value N satisfies

0 <= N < cols.extent(0),

then cols.extent(0) - N of the entries of cols are not valid local column indices. If the returned value is Teuchos::OrdinalTraits<LocalOrdinal>::invalid(), then at least one of the following is true:

  • ! isFillActive ()
  • ! hasColMap ()
  • cols.extent (0) != vals.extent (0)

Definition at line 2401 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceLocalValues ( const LocalOrdinal  localRow,
const Teuchos::ArrayView< const LocalOrdinal > &  cols,
const Teuchos::ArrayView< const Scalar > &  vals 
)

Backwards compatibility version of replaceLocalValues (see above), that takes Teuchos::ArrayView (host pointers) instead of Kokkos::View.

Definition at line 2381 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceLocalValues ( const LocalOrdinal  localRow,
const LocalOrdinal  numEnt,
const Scalar  inputVals[],
const LocalOrdinal  inputCols[] 
)

Epetra compatibility version of replaceLocalValues, that takes raw pointers instead of Kokkos::View.

This version of the method takes the same arguments in the same order as Epetra_CrsMatrix::ReplaceMyValues.

Parameters
localRow[in] local index of the row in which to replace the entries. This row must be owned by the calling process.
numEnt[in] Number of entries to replace; number of valid entries in inputVals and inputCols.
inputVals[in] Values to use for replacing the entries.
inputCols[in] Local indices of the columns in which to replace the entries.
Returns
The number of indices for which values were actually replaced; the number of "correct" indices.

Definition at line 2420 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal, class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoGlobalValuesImpl ( impl_scalar_type  rowVals[],
const crs_graph_type graph,
const RowInfo rowInfo,
const GlobalOrdinal  inds[],
const impl_scalar_type  newVals[],
const LocalOrdinal  numElts,
const bool  atomic = useAtomicUpdatesByDefault 
)
protectedvirtual

Implementation detail of sumIntoGlobalValues.

Template Parameters
InputMemorySpaceKokkos memory space / device in which the input data live. This may differ from the memory space in which the current matrix values (rowVals) live.
ValsMemorySpaceKokkos memory space / device in which the matrix's current values live. This may differ from the memory space in which the input data (inds and newVals) live.
Parameters
rowVals[in/out] On input: Values of the row of the sparse matrix to modify. On output: The modified values.
graph[in] The matrix's graph.
rowInfo[in] Result of getRowInfo on the index of the local row of the matrix to modify.
inds[in] Global column indices of that row to modify.
newVals[in] For each k, increment the value in rowVals corresponding to global column index inds[k] by newVals[k].
Returns
The number of valid input column indices. In case of error other than one or more invalid column indices, this method returns Teuchos::OrdinalTraits<LocalOrdinal>::invalid().

Definition at line 2542 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal, class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoGlobalValues ( const GlobalOrdinal  globalRow,
const Teuchos::ArrayView< const GlobalOrdinal > &  cols,
const Teuchos::ArrayView< const Scalar > &  vals,
const bool  atomic = useAtomicUpdatesByDefault 
)

Sum into one or more sparse matrix entries, using global indices.

This is a local operation; it does not involve communication. However, if you sum into rows not owned by the calling process, it may result in future communication in globalAssemble() (which is called by fillComplete()).

If globalRow is owned by the calling process, then this method performs the sum-into operation right away. Otherwise, if the row is not owned by the calling process, this method defers the sum-into operation until globalAssemble(). That method communicates data for nonowned rows to the processes that own those rows. Then, globalAssemble() does one of the following:

  • It calls insertGlobalValues() for that data if the matrix has a dynamic graph.
  • It calls sumIntoGlobalValues() for that data if the matrix has a static graph. The matrix silently ignores (row,column) pairs that do not exist in the graph.
Parameters
globalRow[in] The global index of the row in which to sum into the matrix entries.
cols[in] One or more column indices.
vals[in] One or more values corresponding to those column indices. vals[k] corresponds to cols[k].
atomic[in] Whether to use atomic updates.
Returns
The number of indices for which values were actually modified; the number of "correct" indices.

This method has the same preconditions and return value meaning as replaceGlobalValues() (which see).

Definition at line 2626 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoGlobalValues ( const GlobalOrdinal  globalRow,
const LocalOrdinal  numEnt,
const Scalar  vals[],
const GlobalOrdinal  cols[],
const bool  atomic = useAtomicUpdatesByDefault 
)

Epetra compatibility version of sumIntoGlobalValues (see above), that takes input as raw pointers instead of Kokkos::View.

Arguments are the same and in the same order as those of Epetra_CrsMatrix::SumIntoGlobalValues, except for atomic, which is as above.

Parameters
globalRow[in] The global index of the row in which to sum into the matrix entries.
numEnt[in] Number of valid entries in vals and cols. This has type LocalOrdinal because we assume that users will never want to insert more column indices in one call than the matrix has columns.
vals[in] numEnt values corresponding to the column indices in cols. That is, vals[k] is the value corresponding to cols[k].
cols[in] numEnt global column indices.
atomic[in] Whether to use atomic updates.
Returns
The number of indices for which values were actually modified; the number of "correct" indices.

Definition at line 2646 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal , class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoLocalValuesImpl ( impl_scalar_type  rowVals[],
const crs_graph_type graph,
const RowInfo rowInfo,
const LocalOrdinal  inds[],
const impl_scalar_type  newVals[],
const LocalOrdinal  numElts,
const bool  atomic = useAtomicUpdatesByDefault 
)
protectedvirtual

Implementation detail of sumIntoLocalValues.

Parameters
rowVals[in/out] On input: Values of the row of the sparse matrix to modify. On output: The modified values.
graph[in] The matrix's graph.
rowInfo[in] Result of graph.getRowInfo on the index of the local row of the matrix to modify.
inds[in] Local column indices of that row to modify.
newVals[in] For each k, increment the value in rowVals corresponding to local column index inds[k] by newVals[k].
atomic[in] Whether to use atomic updates (+=) when incrementing values.

Definition at line 2980 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoLocalValues ( const local_ordinal_type  localRow,
const Kokkos::View< const local_ordinal_type *, Kokkos::AnonymousSpace > &  inputInds,
const Kokkos::View< const impl_scalar_type *, Kokkos::AnonymousSpace > &  inputVals,
const bool  atomic = useAtomicUpdatesByDefault 
)

Sum into one or more sparse matrix entries, using local row and column indices.

For local row index localRow and local column indices cols, perform the update A(localRow, cols[k]) += vals[k]. The row index and column indices must be valid on the calling process, and all matrix entries A(localRow, cols[k]) must already exist. (This method does not change the matrix's structure.) If the row index is valid, any invalid column indices are ignored, but counted in the return value.

This overload of the method takes the column indices and values as Kokkos::View. See below for an overload that takes Teuchos::ArrayView instead.

Template Parameters
LocalIndicesViewTypeKokkos::View specialization that is a 1-D array of LocalOrdinal.
ImplScalarViewTypeKokkos::View specialization that is a 1-D array of impl_scalar_type (usually the same as Scalar, unless Scalar is std::complex<T> for some T, in which case it is Kokkos::complex<T>).
Parameters
localRow[in] Local index of a row. This row must be owned by the calling process.
cols[in] Local indices of the columns whose entries we want to modify.
vals[in] Values corresponding to the above column indices. vals(k) corresponds to cols(k).
atomic[in] Whether to use atomic updates.
Returns
The number of indices for which values were actually modified; the number of "correct" indices.

This method has the same preconditions and return value meaning as replaceLocalValues() (which see).

Definition at line 3082 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoLocalValues ( const LocalOrdinal  localRow,
const Teuchos::ArrayView< const LocalOrdinal > &  cols,
const Teuchos::ArrayView< const Scalar > &  vals,
const bool  atomic = useAtomicUpdatesByDefault 
)

Sum into one or more sparse matrix entries, using local row and column indices.

For local row index localRow and local column indices cols, perform the update A(localRow, cols[k]) += vals[k]. The row index and column indices must be valid on the calling process, and all matrix entries A(localRow, cols[k]) must already exist. (This method does not change the matrix's structure.) If the row index is valid, any invalid column indices are ignored, but counted in the return value.

This overload of the method takes the column indices and values as Teuchos::ArrayView. See above for an overload that takes Kokkos::View instead.

Parameters
localRow[in] Local index of a row. This row must be owned by the calling process.
cols[in] Local indices of the columns whose entries we want to modify.
vals[in] Values corresponding to the above column indices. vals[k] corresponds to cols[k].
atomic[in] Whether to use atomic updates.
Returns
The number of indices for which values were actually modified; the number of "correct" indices.

This method has the same preconditions and return value meaning as replaceLocalValues() (which see).

Definition at line 3062 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sumIntoLocalValues ( const LocalOrdinal  localRow,
const LocalOrdinal  numEnt,
const Scalar  vals[],
const LocalOrdinal  cols[],
const bool  atomic = useAtomicUpdatesByDefault 
)

Epetra compatibility version of sumIntoLocalValues (see above) that takes raw pointers instead of Kokkos::View.

Arguments are the same and in the same order as Epetra_CrsMatrix::SumIntoMyValues, except for the atomic last argument, which is as above.

Parameters
localRow[in] The local index of the row in which to sum into the matrix entries.
numEnt[in] Number of valid entries in vals and cols. This has type LocalOrdinal because we assume that users will never want to insert more column indices in one call than the matrix has columns.
vals[in] numEnt values corresponding to the column indices in cols. That is, vals[k] is the value corresponding to cols[k].
cols[in] numEnt local column indices.
atomic[in] Whether to use atomic updates.
Returns
The number of indices for which values were actually modified; the number of "correct" indices.

Definition at line 3102 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
template<class LocalIndicesViewType , class ImplScalarViewType , class BinaryFunction >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::transformLocalValues ( const LocalOrdinal  lclRow,
const typename UnmanagedView< LocalIndicesViewType >::type &  inputInds,
const typename UnmanagedView< ImplScalarViewType >::type &  inputVals,
BinaryFunction  f,
const bool  atomic = useAtomicUpdatesByDefault 
)
inline

Transform CrsMatrix entries in place, using local indices to select the entries in the row to transform.

For every entry $A(i,j)$ to transform, if $v_{ij}$ is the corresponding entry of the inputVals array, then we apply the binary function f to $A(i,j)$ as follows:

\[ A(i,j) := f(A(i,j), v_{ij}). \]

For example, BinaryFunction = std::plus<impl_scalar_type> does the same thing as sumIntoLocalValues, and BinaryFunction = project2nd<impl_scalar_type,impl_scalar_type> does the same thing as replaceLocalValues. (It is generally more efficient to call sumIntoLocalValues resp. replaceLocalValues than to do this.)

This overload of the method takes the column indices and values as Kokkos::View. See below for an overload that takes Teuchos::ArrayView instead.

Template Parameters
LocalIndicesViewTypeKokkos::View specialization that is a 1-D array of LocalOrdinal.
ImplScalarViewTypeKokkos::View specialization that is a 1-D array of impl_scalar_type (usually the same as Scalar, unless Scalar is std::complex<T> for some T, in which case it is Kokkos::complex<T>).
BinaryFunctionThe type of the binary function f to use for updating the sparse matrix's value(s). This should be convertible to std::function<impl_scalar_type (const impl_scalar_type&, const impl_scalar_type&)>.
Parameters
lclRow[in] (Local) index of the row to modify. This row must</t> be owned by the calling process. (This is a stricter requirement than for sumIntoGlobalValues.)
inputInds[in] (Local) indices in the row to modify. Indices not in the row on the calling process, and their corresponding values, will be ignored.
inputVals[in] Values to use for modification.
f[in] The binary function to use for updating the sparse matrix's value. It takes two impl_scalar_type values and returns impl_scalar_type.
atomic[in] Whether to use atomic updates.

Definition at line 1709 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
template<class BinaryFunction , class InputMemorySpace >
LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::transformGlobalValues ( const GlobalOrdinal  gblRow,
const Kokkos::View< const GlobalOrdinal *, InputMemorySpace, Kokkos::MemoryUnmanaged > &  inputInds,
const Kokkos::View< const impl_scalar_type *, InputMemorySpace, Kokkos::MemoryUnmanaged > &  inputVals,
BinaryFunction  f,
const bool  atomic = useAtomicUpdatesByDefault 
)
inline

Transform CrsMatrix entries in place, using global indices to select the entries in the row to transform.

For every entry $A(i,j)$ to transform, if $v_{ij}$ is the corresponding entry of the inputVals array, then we apply the binary function f to $A(i,j)$ as follows:

\[ A(i,j) := f(A(i,j), v_{ij}). \]

For example, BinaryFunction = std::plus<impl_scalar_type> does the same thing as sumIntoLocalValues, and BinaryFunction = project2nd<impl_scalar_type,impl_scalar_type> does the same thing as replaceLocalValues. (It is generally more efficient to call sumIntoLocalValues resp. replaceLocalValues than to do this.)

Template Parameters
BinaryFunctionThe type of the binary function f to use for updating the sparse matrix's value(s). This should be convertible to std::function<impl_scalar_type (const impl_scalar_type&, const impl_scalar_type&)>.
InputMemorySpaceKokkos memory space / device in which the input data live. This may differ from the memory space in which the current matrix's row's values live.
Parameters
gblRow[in] (Global) index of the row to modify. This row must</t> be owned by the calling process. (This is a stricter requirement than for sumIntoGlobalValues.)
inputInds[in] (Global) indices in the row to modify. Indices not in the row on the calling process, and their corresponding values, will be ignored.
inputVals[in] Values to use for modification.
f[in] The binary function to use for updating the sparse matrix's value. It takes two impl_scalar_type values and returns impl_scalar_type.
atomic[in] Whether to use atomic updates.

This method works whether indices are local or global. However, it will cost more if indices are local, since it will have to convert the input global indices to local indices in that case.

Definition at line 1799 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setAllToScalar ( const Scalar &  alpha)

Set all matrix entries equal to alpha.

Definition at line 3427 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scale ( const Scalar &  alpha)

Scale the matrix's values: this := alpha*this.

Definition at line 3406 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setAllValues ( const typename local_graph_device_type::row_map_type &  ptr,
const typename local_graph_device_type::entries_type::non_const_type &  ind,
const typename local_matrix_device_type::values_type &  val 
)

Set the local matrix using three (compressed sparse row) arrays.

Precondition
ind is sorted within each row
hasColMap() == true
getGraph() != Teuchos::null
No insert/sum routines have been called
Warning
This is for EXPERT USE ONLY. We make NO PROMISES of backwards compatibility.

This method behaves like the CrsMatrix constructor that takes a const CrsGraph. It fixes the matrix's graph, but does not call fillComplete on the matrix. The graph might not necessarily be fill complete, but it must have a local graph.

The input arguments might be used directly (shallow copy), or they might be (deep) copied.

Parameters
ptr[in] Array of row offsets.
ind[in] Array of (local) column indices.
val[in/out] Array of values. This is in/out because the matrix reserves the right to take this argument by shallow copy. Any method that changes the matrix's values may then change this.

Definition at line 3451 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setAllValues ( const local_matrix_device_type localMatrix)

Set the local matrix using an existing local matrix.

Precondition
column indices are sorted within each row
hasColMap() == true
getGraph() != Teuchos::null
No insert/sum routines have been called
Warning
This is for EXPERT USE ONLY. We make NO PROMISES of backwards compatibility.

This method simply calls the method setAllValues that accepts three compressed sparse row arrays.

The input argument might be used directly (shallow copy), or it might be (deep) copied.

Parameters
ptr[in/out] Kokkos sparse local matrix This is in/out because the matrix reserves the right to take this argument by shallow copy. Any method that changes the matrix's values may then change this.

Definition at line 3499 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setAllValues ( const Teuchos::ArrayRCP< size_t > &  ptr,
const Teuchos::ArrayRCP< LocalOrdinal > &  ind,
const Teuchos::ArrayRCP< Scalar > &  val 
)

Set the local matrix using three (compressed sparse row) arrays.

Precondition
ind is sorted within each row
hasColMap() == true
getGraph() != Teuchos::null
No insert/sum routines have been called
Warning
This is for EXPERT USE ONLY. We make NO PROMISES of backwards compatibility.

This method behaves like the CrsMatrix constructor that takes a const CrsGraph. It fixes the matrix's graph, but does not call fillComplete on the matrix. The graph might not necessarily be fill complete, but it must have a local graph.

The input arguments might be used directly (shallow copy), or they might be (deep) copied.

Parameters
ptr[in] Array of row offsets.
ind[in] Array of (local) column indices.
val[in/out] Array of values. This is in/out because the matrix reserves the right to take this argument by shallow copy. Any method that changes the matrix's values may then change this.

Definition at line 3517 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
row_ptrs_host_view_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowPtrsHost ( ) const
inline

Get a host view of the CRS packed row pointers.

Definition at line 1910 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
row_ptrs_device_view_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowPtrsDevice ( ) const
inline

Get a device view of the CRS packed row pointers.

Definition at line 1914 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
local_inds_host_view_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalIndicesHost ( ) const
inline

Get a host view of the CRS packed column indicies.

Definition at line 1918 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
local_inds_device_view_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalIndicesDevice ( ) const
inline

Get a device_view of the CRS packed column indicies.

Definition at line 1922 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::globalAssemble ( )

Communicate nonlocal contributions to other processes.

Users do not normally need to call this method. fillComplete always calls this method, unless you specifically tell fillComplete to do otherwise by setting its "No Nonlocal Changes" parameter to true. Thus, it suffices to call fillComplete.

Methods like insertGlobalValues and sumIntoGlobalValues let you add or modify entries in rows that are not owned by the calling process. These entries are called "nonlocal contributions." The methods that allow nonlocal contributions store the entries on the calling process, until globalAssemble is called. globalAssemble sends these nonlocal contributions to the process(es) that own them, where they then become part of the matrix.

This method only does global assembly if there are nonlocal entries on at least one process. It does an all-reduce to find that out. If not, it returns early, without doing any more communication or work.

If you previously inserted into a row which is not owned by any process in the row Map, the behavior of this method is undefined. It may detect the invalid row indices and throw an exception, or it may silently drop the entries inserted into invalid rows. Behavior may vary, depending on whether Tpetra was built with debug checking enabled.

Definition at line 4056 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::resumeFill ( const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null)

Resume operations that may change the values or structure of the matrix.

This method must be called as a collective operation.

Calling fillComplete "freezes" both the values and the structure of the matrix. If you want to modify the matrix again, you must first call resumeFill. You then may not call resumeFill again on that matrix until you first call fillComplete. You may make sequences of fillComplete, resumeFill calls as many times as you wish.

Postcondition
isFillActive() && ! isFillComplete()

Definition at line 4305 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fillComplete ( const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Tell the matrix that you are done changing its structure or values, and that you are ready to do computational kernels (e.g., sparse matrix-vector multiply) with it.

This tells the graph to optimize its data structures for computational kernels, and to prepare (MPI) communication patterns.

Off-process indices are distributed (via globalAssemble()), indices are sorted, redundant indices are fused, and global indices are transformed to local indices.

Warning
The domain Map and row Map arguments to this method MUST be one to one! If you have Maps that are not one to one, and you do not know how to make a Map that covers the same global indices but is one to one, then you may call Tpetra::createOneToOne() (see Map's header file) to make a one-to-one version of your Map.
Precondition
isFillActive() && ! isFillComplete()
Postcondition
! isFillActive() && isFillComplete()
Parameters
domainMap[in] The matrix's domain Map. MUST be one to one!
rangeMap[in] The matrix's range Map. MUST be one to one! May be, but need not be, the same as the domain Map.
params[in/out] List of parameters controlling this method's behavior. See below for valid parameters.

List of valid parameters in params:

  • "No Nonlocal Changes" (bool): Default is false. If true, the caller promises that no modifications to nonowned rows have happened on any process since the last call to fillComplete. This saves a global all-reduce to check whether any process did a nonlocal insert. Nonlocal changes include any sumIntoGlobalValues or insertGlobalValues call with a row index that is not in the row Map of the calling process.

  • "Sort column Map ghost GIDs" (bool): Default is true. makeColMap() (which fillComplete may call) always groups remote GIDs by process rank, so that all remote GIDs with the same owning rank occur contiguously. By default, it always sorts remote GIDs in increasing order within those groups. This behavior differs from Epetra, which does not sort remote GIDs with the same owning process. If you don't want to sort (for compatibility with Epetra), set this parameter to false. This parameter only takes effect if the matrix owns the graph. This is an expert mode parameter ONLY. We make no promises about backwards compatibility of this parameter. It may change or disappear at any time.

Definition at line 4349 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fillComplete ( const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null)

Tell the matrix that you are done changing its structure or values, and that you are ready to do computational kernels (e.g., sparse matrix-vector multiply) with it. Set default domain and range Maps.

See above three-argument version of fillComplete for full documentation. If the matrix does not yet have domain and range Maps (i.e., if fillComplete has not yet been called on this matrix at least once), then this method uses the matrix's row Map (result of this->getRowMap()) as both the domain Map and the range Map. Otherwise, this method uses the matrix's existing domain and range Maps.

Warning
It is only valid to call this overload of fillComplete if the row Map is one to one! If the row Map is NOT one to one, you must call the above three-argument version of fillComplete, and supply one-to-one domain and range Maps. If you have Maps that are not one to one, and you do not know how to make a Map that covers the same global indices but is one to one, then you may call Tpetra::createOneToOne() (see Map's header file) to make a one-to-one version of your Map.
Parameters
params[in/out] List of parameters controlling this method's behavior. See documentation of the three-argument version of fillComplete (above) for valid parameters.

Definition at line 4323 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::expertStaticFillComplete ( const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const Teuchos::RCP< const import_type > &  importer = Teuchos::null,
const Teuchos::RCP< const export_type > &  exporter = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Perform a fillComplete on a matrix that already has data.

The matrix must already have filled local 1-D storage (k_clInds1D_ and k_rowPtrs_ for the graph, and k_values1D_ in the matrix). If the matrix has been constructed in any other way, this method will throw an exception. This routine is needed to support other Trilinos packages and should not be called by ordinary users.

Warning
This method is intended for expert developer use only, and should never be called by user code.
Parameters
domainMap[in] The matrix's domain Map. MUST be one to one!
rangeMap[in] The matrix's range Map. MUST be one to one! May be, but need not be, the same as the domain Map.
importer[in] Import from the matrix's domain Map to its column Map. If no Import is necessary (i.e., if the domain and column Maps are the same, in the sense of Tpetra::Map::isSameAs), then this may be Teuchos::null.
exporter[in] Export from the matrix's row Map to its range Map. If no Export is necessary (i.e., if the row and range Maps are the same, in the sense of Tpetra::Map::isSameAs), then this may be Teuchos::null.
params[in/out] List of parameters controlling this method's behavior.

Definition at line 4548 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceColMap ( const Teuchos::RCP< const map_type > &  newColMap)

Replace the matrix's column Map with the given Map.

Parameters
newColMap[in] New column Map. Must be nonnull. Within Tpetra, there are no particular restrictions on the column map. However, if this graph will be used in Xpetra, Ifpack2, or MueLu, the column map's list of global indices must follow "Aztec ordering": locally owned GIDs (same order as in domain map), followed by remote GIDs (in order of owning proc, and sorted within each proc).

It is strongly recommended to use Tpetra::Details::makeColMap() to create the column map. makeColMap() follows Aztec ordering by default.

Precondition
The matrix must have no entries inserted yet, on any process in the row Map's communicator.
The matrix must not have been created with a constant (a.k.a. "static") CrsGraph.

Definition at line 3914 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::reindexColumns ( crs_graph_type *const  graph,
const Teuchos::RCP< const map_type > &  newColMap,
const Teuchos::RCP< const import_type > &  newImport = Teuchos::null,
const bool  sortEachRow = true 
)

Reindex the column indices in place, and replace the column Map. Optionally, replace the Import object as well.

Precondition
The matrix is not fill complete: ! this->isFillComplete() .
Either the input graph is NULL, or it is not fill complete: graph == NULL || ! graph->isFillComplete().
On every calling process, every index owned by the current column Map must also be owned by the new column Map.
If the new Import object is provided, the new Import object's source Map must be the same as the current domain Map, and the new Import's target Map must be the same as the new column Map.
Parameters
graph[in] The matrix's graph. If you don't provide this (i.e., if graph == NULL), then the matrix must own its graph, which will be modified in place. (That is, you must not have created the matrix with a constant graph.) If you do provide this, then the method will assume that it is the same graph as the matrix's graph, and the provided graph will be modified in place.
newColMap[in] New column Map. Must be nonnull.
newImport[in] New Import object. Optional; computed if not provided or if null. Computing an Import is expensive, so it is worth providing this if you can.
sortEachRow[in] If true, sort the indices (and their corresponding values) in each row after reindexing.

Why would you want to use this method? Well, for example, you might need to use an Ifpack2 preconditioner that only accepts a matrix with a certain kind of column Map. Your matrix has the wrong kind of column Map, but you know how to compute the right kind of column Map. You might also know an efficient way to compute an Import object from the current domain Map to the new column Map. (For an instance of the latter, see the Details::makeOptimizedColMapAndImport function in Tpetra_Details_makeOptimizedColMap.hpp.)

Suppose that you created this CrsMatrix with a constant graph; that is, that you called the CrsMatrix constructor that takes a CrsGraph as input:

RCP<CrsGraph<> > G (new CrsGraph<> (rowMap, origColMap, ...));
// ... fill G ...
G->fillComplete (domMap, ranMap);
CrsMatrix<> A (G);
// ... fill A ...

Now suppose that you want to give A to a preconditioner that can't handle a matrix with an arbitrary column Map (in the example above, origColMap). You first must create a new suitable column Map newColMap, and optionally a new Import object newImport from the matrix's current domain Map to the new column Map. Then, call this method, passing in G (which must not be fill complete) while the matrix is not fill complete. Be sure to save the graph's original Import object; you'll need that later.

RCP<const CrsGraph<>::import_type> origImport = G->getImporter ();
G->resumeFill ();
A.reindexColumns (G.getRawPtr (), newColMap, newImport);
G.fillComplete (domMap, ranMap);
A.fillComplete (domMap, ranMap);

Now you may give the matrix A to the preconditioner in question. After doing so, and after you solve the linear system using the preconditioner, you might want to put the matrix back like it originally was. You can do that, too!

A.resumeFill ();
G->resumeFill ();
A.reindexColumns (G.getRawPtr (), origColMap, origImport);
G->fillComplete (domMap, ranMap);
A->fillComplete (domMap, ranMap);

Definition at line 3932 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceDomainMap ( const Teuchos::RCP< const map_type > &  newDomainMap)

Replace the current domain Map with the given objects.

The matrix's Import object will be recomputed if needed.

Parameters
newDomainMap[in] New domain Map. Must be nonnull.
Precondition
The matrix must be fill complete: isFillComplete() == true.

Definition at line 3968 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceDomainMapAndImporter ( const Teuchos::RCP< const map_type > &  newDomainMap,
Teuchos::RCP< const import_type > &  newImporter 
)

Replace the current domain Map and Import with the given objects.

Parameters
newDomainMap[in] New domain Map. Must be nonnull.
newImporter[in] Optional Import object. If null, the new Domain Map must equal the matrix's Column Map
Precondition
The matrix must be fill complete: isFillComplete() == true.
If the Import is provided, its target Map must be the same as the column Map of the matrix.
If the Import is provided, its source Map must be the same as the provided new domain Map.
If the Import is not provided, the new Domain Map must be the same as the matrix's Column Map.

Definition at line 3983 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceRangeMap ( const Teuchos::RCP< const map_type > &  newRangeMap)

Replace the current range Map with the given objects.

The matrix's Export object will be recomputed if needed.

Parameters
newRangeMap[in] New Range Map. Must be nonnull.
Precondition
The matrix must be fill complete: isFillComplete() == true.

Definition at line 3999 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::replaceRangeMapAndExporter ( const Teuchos::RCP< const map_type > &  newRangeMap,
Teuchos::RCP< const export_type > &  newExporter 
)

Replace the current Range Map and Export with the given objects.

Parameters
newRangeMap[in] New domain Map. Must be nonnull.
newExporter[in] Optional Export object. If null, the new Range Map must equal the matrix's Row Map
Precondition
The matrix must be fill complete: isFillComplete() == true.
If the Export is provided, its target Map must be the same as the new Range Map of the matrix.
If the Export is provided, its source Map must be the same as the Row Map of the matrix
If the Export is not provided, the new Range Map must be the same as the matrix's Row Map.

Definition at line 4014 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::removeEmptyProcessesInPlace ( const Teuchos::RCP< const map_type > &  newMap)
overridevirtual

Remove processes owning zero rows from the Maps and their communicator.

Warning
This method is ONLY for use by experts. We highly recommend using the nonmember function of the same name defined in Tpetra_DistObject_decl.hpp.
We make NO promises of backwards compatibility. This method may change or disappear at any time.
Parameters
newMap[in] This must be the result of calling the removeEmptyProcesses() method on the row Map. If it is not, this method's behavior is undefined. This pointer will be null on excluded processes.

Definition at line 7378 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Teuchos::Comm< int > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getComm ( ) const
overridevirtual

The communicator over which the matrix is distributed.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 809 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRowMap ( ) const
overridevirtual

The Map that describes the row distribution in this matrix.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 937 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getColMap ( ) const
overridevirtual

The Map that describes the column distribution in this matrix.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 944 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const RowGraph< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGraph ( ) const
overridevirtual

This matrix's graph, as a RowGraph.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 965 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getCrsGraph ( ) const

This matrix's graph, as a CrsGraph.

Definition at line 975 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_device_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalMatrixDevice ( ) const

The local sparse matrix.

Warning
It is only valid to call this method under certain circumstances. In particular, either the CrsMatrix must have been created with a local_matrix_type object, or fillComplete must have been called on this CrsMatrix at least once. This method will do no error checking, so you are responsible for knowing when it is safe to call this method.

Definition at line 1011 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
std::shared_ptr< typename CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_multiply_op_type > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalMultiplyOperator ( ) const

The local sparse matrix operator (a wrapper of getLocalMatrixDevice() that supports local matrix-vector multiply)

Warning
It is only valid to call this method if this->isFillComplete().

Definition at line 1035 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
global_size_t Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumRows ( ) const
overridevirtual

Number of global elements in the row map of this matrix.

This is <it>not</it> the number of rows in the matrix as a mathematical object. This method returns the global sum of the number of local elements in the row map on each processor, which is the row map's getGlobalNumElements(). Since the row map is not one-to-one in general, that global sum could be different than the number of rows in the matrix. If you want the number of rows in the matrix, ask the range map for its global number of elements, using the following code: global_size_t globalNumRows = getRangeMap()->getGlobalNumElements(); This method retains the behavior of Epetra, which also asks the row map for the global number of rows, rather than asking the range map.

Warning
Undefined if isFillActive().

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 872 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
global_size_t Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumCols ( ) const
overridevirtual

The number of global columns in the matrix.

This equals the number of entries in the matrix's domain Map.

Warning
Undefined if isFillActive().

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 879 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalNumRows ( ) const
overridevirtual

The number of matrix rows owned by the calling process.

Note that the sum of all the return values over all processes in the row Map's communicator does not necessarily equal the global number of rows in the matrix, if the row Map is overlapping.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 886 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalNumCols ( ) const
overridevirtual

The number of columns connected to the locally owned rows of this matrix.

Throws std::runtime_error if ! hasColMap ().

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 894 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
GlobalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getIndexBase ( ) const
overridevirtual

The index base for global indices for this matrix.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 930 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
global_size_t Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumEntries ( ) const
overridevirtual

The global number of entries in this matrix.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 858 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalNumEntries ( ) const
overridevirtual

The local number of entries in this matrix.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 865 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal, class Node >
size_t Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInGlobalRow ( GlobalOrdinal  globalRow) const
overridevirtual

Number of entries in the sparse matrix in the given global row, on the calling (MPI) process.

Returns
OrdinalTraits<size_t>::invalid()if the specified global row index is invalid on the calling process, else the number of entries in the given row.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 902 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInLocalRow ( local_ordinal_type  localRow) const
overridevirtual

Number of entries in the sparse matrix in the given local row, on the calling (MPI) process.

Returns
OrdinalTraits<size_t>::invalid()if the specified local row index is invalid on the calling process, else the number of entries in the given row.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 909 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalMaxNumRowEntries ( ) const
overridevirtual

Maximum number of entries in any row of the matrix, over all processes in the matrix's communicator.

Precondition
! isFillActive()

This method only uses the matrix's graph. Explicitly stored zeros count as "entries."

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 916 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalMaxNumRowEntries ( ) const
overridevirtual

Maximum number of entries in any row of the matrix, on this process.

Precondition
! isFillActive()

This method only uses the matrix's graph. Explicitly stored zeros count as "entries."

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 923 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual LocalOrdinal Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getBlockSize ( ) const
inlineoverridevirtual

The number of degrees of freedom per mesh point.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 2418 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::hasColMap ( ) const
overridevirtual

Whether the matrix has a well-defined column Map.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 851 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isLocallyIndexed ( ) const
overridevirtual

Whether the matrix is locally indexed on the calling process.

The matrix is locally indexed on the calling process if and only if all of the following hold:

  1. The matrix is not empty on the calling process
  2. The matrix has a column Map

The following is always true:

(! locallyIndexed() && ! globallyIndexed()) || (locallyIndexed() || globallyIndexed());

That is, a matrix may be neither locally nor globally indexed, but it can never be both. Furthermore a matrix that is not fill complete, might have some processes that are neither locally nor globally indexed, and some processes that are globally indexed. The processes that are neither do not have any entries.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 837 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isGloballyIndexed ( ) const
overridevirtual

Whether the matrix is globally indexed on the calling process.

The matrix is globally indexed on the calling process if and only if all of the following hold:

  1. The matrix is not empty on the calling process
  2. The matrix does not yet have a column Map

The following is always true:

(! locallyIndexed() && ! globallyIndexed()) ||
(locallyIndexed() || globallyIndexed());

That is, a matrix may be neither locally nor globally indexed, but it can never be both. Furthermore a matrix that is not fill complete, might have some processes that are neither locally nor globally indexed, and some processes that are globally indexed. The processes that are neither do not have any entries.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 844 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isFillComplete ( ) const
overridevirtual

Whether the matrix is fill complete.

A matrix is fill complete (or "in compute mode") when fillComplete() has been called without an intervening call to resumeFill(). A matrix must be fill complete in order to call computational kernels like sparse matrix-vector multiply and sparse triangular solve. A matrix must be not fill complete ("in edit mode") in order to call methods that insert, modify, or remove entries.

The following are always true:

A matrix starts out (after its constructor returns) as not fill complete. It becomes fill complete after fillComplete() returns, and becomes not fill complete again if resumeFill() is called. Some methods like clone() and some of the "nonmember constructors" (like importAndFillComplete() and exportAndFillComplete()) may return a fill-complete matrix.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 816 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isFillActive ( ) const

Whether the matrix is not fill complete.

A matrix is fill complete (or "in compute mode") when fillComplete() has been called without an intervening call to resumeFill(). A matrix must be fill complete in order to call computational kernels like sparse matrix-vector multiply and sparse triangular solve. A matrix must be not fill complete ("in edit mode") in order to call methods that insert, modify, or remove entries.

The following are always true:

A matrix starts out (after its constructor returns) as not fill complete. It becomes fill complete after fillComplete() returns, and becomes not fill complete again if resumeFill() is called. Some methods like clone() and some of the "nonmember constructors" (like importAndFillComplete() and exportAndFillComplete()) may return a fill-complete matrix.

Definition at line 823 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isStorageOptimized ( ) const

Returns true if storage has been optimized.

Optimized storage means that the allocation of each row is equal to the number of entries. The effect is that a pass through the matrix, i.e., during a mat-vec, requires minimal memory traffic. One limitation of optimized storage is that no new indices can be added to the matrix.

Definition at line 830 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isStaticGraph ( ) const

Indicates that the graph is static, so that new entries cannot be added to this matrix.

Definition at line 1072 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getFrobeniusNorm ( ) const
overridevirtual

Compute and return the Frobenius norm of the matrix.

The Frobenius norm of the matrix is defined as [ |A|_F = {{i,j} |A(i,j)|^2}. ].

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3862 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::supportsRowViews ( ) const
overridevirtual

Return true if getLocalRowView() and getGlobalRowView() are valid for this object.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 1086 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal, class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowCopy ( GlobalOrdinal  GlobalRow,
nonconst_global_inds_host_view_type &  Indices,
nonconst_values_host_view_type &  Values,
size_t &  NumEntries 
) const
override

Fill given arrays with a deep copy of the locally owned entries of the matrix in a given row, using global column indices.

Parameters
GlobalRow[in] Global index of the row for which to return entries.
Indices[out] Global column indices corresponding to values.
Values[out] Matrix values.
NumEntries[out] Number of entries.
Note
To Tpetra developers: Discussion of whether to use Scalar or impl_scalar_type for output array of matrix values.

If Scalar differs from impl_scalar_type, as for example with std::complex<T> and Kokkos::complex<T>, we must choose which type to use. We must make the same choice as RowMatrix does, else CrsMatrix won't compile, because it won't implement a pure virtual method. We choose Scalar, for the following reasons. First, Scalar is the user's preferred type, and impl_scalar_type an implementation detail that makes Tpetra work with Kokkos. Second, Tpetra's public interface provides a host-only interface, which eliminates some reasons for requiring implementation-specific types like Kokkos::complex.

We do eventually want to put Tpetra methods in Kokkos kernels, but we only need to put them in host kernels, since Tpetra is a host-only interface. Users can still manually handle conversion from Scalar to impl_scalar_type for reductions.

The right thing to do would be to rewrite RowMatrix so that getGlobalRowCopy is NOT inherited, but is implemented by a pure virtual "hook" getGlobalRowCopyImpl. The latter takes raw pointers. That would give us the freedom to overload getGlobalRowCopy, which one normally can't do with virtual methods. It would make sense for one getGlobalRowCopyImpl method to implement both Teuchos::ArrayView and Kokos::View versions of getGlobalRowCopy.

Note: A std::runtime_error exception is thrown if either Indices or Values is not large enough to hold the data associated with row GlobalRow. If row GlobalRow is not owned by the calling process, then Indices and Values are unchanged and NumIndices is returned as Teuchos::OrdinalTraits<size_t>::invalid().

Definition at line 3243 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowCopy ( LocalOrdinal  LocalRow,
nonconst_local_inds_host_view_type &  Indices,
nonconst_values_host_view_type &  Values,
size_t &  NumEntries 
) const
override

Fill given arrays with a deep copy of the locally owned entries of the matrix in a given row, using local column indices.

Parameters
LocalRow[in] Local index of the row for which to return entries.
Indices[out] Local column indices corresponding to values.
Values[out] Matrix values.
NumEntries[out] Number of entries returned.

Note: A std::runtime_error exception is thrown if either colInds or vals is not large enough to hold the data associated with row localRow. If row localRow is not owned by the calling process, then colInds and vals are unchanged and numEntries is returned as Teuchos::OrdinalTraits<size_t>::invalid().

Definition at line 3189 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal, class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowView ( GlobalOrdinal  GlobalRow,
global_inds_host_view_type &  indices,
values_host_view_type &  values 
) const
override

Get a constant, nonpersisting view of a row of this matrix, using global row and column indices.

Parameters
GlobalRow[in] Global index of the row to view.
indices[out] On output: view of the global column indices in the row.
values[out] On output: view of the values in the row.
Precondition
isLocallyIndexed () == false
Postcondition
indices.size () == this->getNumEntriesInGlobalRow (GlobalRow)

If GlobalRow is not a valid global row index on the calling process, then indices is set to null.

Definition at line 3349 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowView ( LocalOrdinal  LocalRow,
local_inds_host_view_type &  indices,
values_host_view_type &  values 
) const
override

Get a constant view of a row of this matrix, using local row and column indices.

Parameters
LocalRow[in] Local index of the row to view.
indices[out] On output: view of the local column indices in the row.
values[out] On output: view of the values in the row.
Precondition
isGloballyIndexed () == false
Postcondition
indices.size () == this->getNumEntriesInLocalRow (LocalRow)

If LocalRow is not a valid local row index on the calling process, then indices is set to null.

Definition at line 3290 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalDiagCopy ( Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  diag) const
overridevirtual

Get a constant, nonpersisting view of a row of this matrix, using local row and column indices, with raw pointers.

This overload exists only if Scalar differs from impl_scalar_type. In that case, this overload takes a Scalar pointer. Get a copy of the diagonal entries of the matrix.

This method returns a Vector with the same Map as this matrix's row Map. On each process, it contains the diagonal entries owned by the calling process. If the matrix has an empty row, the diagonal entry contains a zero.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3598 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalDiagOffsets ( Teuchos::ArrayRCP< size_t > &  offsets) const

Get offsets of the diagonal entries in the matrix.

Warning
This method is DEPRECATED. Call CrsGraph::getLocalDiagOffsets, in particular the overload that returns the offsets as a Kokkos::View.
This method is only for expert users.
We make no promises about backwards compatibility for this method. It may disappear or change at any time.
This method must be called collectively. We reserve the right to do extra checking in a debug build that will require collectives.
Precondition
The matrix must be locally indexed (which means that it has a column Map).
All diagonal entries of the matrix's graph must be populated on this process. Results are undefined otherwise.
Postcondition
offsets.size() == getLocalNumRows()

This method creates an array of offsets of the local diagonal entries in the matrix. This array is suitable for use in the two-argument version of getLocalDiagCopy(). However, its contents are not defined in any other context. For example, you should not rely on offsets[i] being the index of the diagonal entry in the views returned by getLocalRowView(). This may be the case, but it need not be. (For example, we may choose to optimize the lookups down to the optimized storage level, in which case the offsets will be computed with respect to the underlying storage format, rather than with respect to the views.)

Calling any of the following invalidates the output array:

If the matrix has a const ("static") graph, and if that graph is fill complete, then the offsets array remains valid through calls to fillComplete() and resumeFill(). "Invalidates" means that you must call this method again to recompute the offsets.

Definition at line 3556 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalDiagCopy ( Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  diag,
const Kokkos::View< const size_t *, device_type, Kokkos::MemoryUnmanaged > &  offsets 
) const

Variant of getLocalDiagCopy() that uses precomputed offsets.

Warning
This method is only for expert users.
We make no promises about backwards compatibility for this method. It may disappear or change at any time.

This method uses the offsets of the diagonal entries, as precomputed by the Kokkos::View overload of getLocalDiagOffsets(), to speed up copying the diagonal of the matrix. The offsets must be recomputed if any of the following methods are called:

If the matrix has a const ("static") graph, and if that graph is fill complete, then the offsets array remains valid through calls to fillComplete() and resumeFill().

Definition at line 3657 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalDiagCopy ( Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  diag,
const Teuchos::ArrayView< const size_t > &  offsets 
) const

Variant of getLocalDiagCopy() that uses precomputed offsets.

Warning
This overload of the method is DEPRECATED. Call the overload above that returns the offsets as a Kokkos::View.
This method is only for expert users.
We make no promises about backwards compatibility for this method. It may disappear or change at any time.

This method uses the offsets of the diagonal entries, as precomputed by getLocalDiagOffsets(), to speed up copying the diagonal of the matrix. The offsets must be recomputed if any of the following methods are called:

If the matrix has a const ("static") graph, and if that graph is fill complete, then the offsets array remains valid through calls to fillComplete() and resumeFill().

Definition at line 3694 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::leftScale ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  x)
overridevirtual

Scale the matrix on the left with the given Vector.

On return, for all entries i,j in the matrix, $A(i,j) = x(i)*A(i,j)$.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3753 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::rightScale ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  x)
overridevirtual

Scale the matrix on the right with the given Vector.

On return, for all entries i,j in the matrix, $A(i,j) = x(j)*A(i,j)$.

Implements Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 3808 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::localApply ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  X,
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  Y,
const Teuchos::ETransp  mode = Teuchos::NO_TRANS,
const Scalar &  alpha = Teuchos::ScalarTraits<Scalar>::one (),
const Scalar &  beta = Teuchos::ScalarTraits<Scalar>::zero () 
) const

Compute the local part of a sparse matrix-(Multi)Vector multiply.

Most Tpetra users want the apply() method (which see), not this method.

This method computes Y := beta*Y + alpha*Op(A)*X, where Op(A) is either $A$, $A^T$ (the transpose), or $A^H$ (the conjugate transpose).

"The local part" means that this method does no communication between processes, even if this is necessary for correctness of the matrix-vector multiply. Use the apply() method if you want to compute mathematical sparse matrix-vector multiply.

This method is mainly of use to Tpetra developers, though some users may find it helpful if they plan to reuse the result of doing an Import on the input MultiVector for several sparse matrix-vector multiplies with matrices that have the same column Map.

If you want to do global (not just local) sparse matrix-vector multiplies with MultiVectors that have different Scalar type than the CrsMatrix, use CrsMatrixMultiplyOp. If you want to do local sparse matrix-vector multiplies with MultiVectors that have different Scalar type than the CrsMatrix (i.e., if you're wondering where the templated localMultiply method went), use LocalCrsMatrixOperator.

The Map of X and mode must satisfy the following:

mode == Teuchos::NO_TRANS &&
X.getMap ()->isSameAs(* (this->getColMap ())) ||
mode != Teuchos::NO_TRANS &&
X.getMap ()->isSameAs(* (this->getRowMap ()));

The Map of Y and mode must satisfy the following:

mode == Teuchos::NO_TRANS &&
Y.getMap ()->isSameAs(* (this->getRowMap ())) ||
mode != Teuchos::NO_TRANS &&
Y.getMap ()->isSameAs(* (this->getColMap ()));

We say that X is "post-Imported," and that Y is "pre-Exported."

Both X and Y must have constant stride, and they may not alias one another (that is, occupy overlapping space in memory). We may not necessarily check for aliasing, and if we do, we will only do this in a debug build. Aliasing X and Y may cause nondeterministically incorrect results.

If beta == 0, this operation will enjoy overwrite semantics: Y's entries will be ignored, and Y will be overwritten with the result of the multiplication, even if it contains Inf or NaN floating-point entries. Likewise, if alpha == 0, this operation will ignore A and X, even if they contain Inf or NaN floating-point entries.

Definition at line 5009 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
template<class T >
Teuchos::RCP< CrsMatrix< T, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::convert ( ) const

Return another CrsMatrix with the same entries, but converted to a different Scalar type T.

Definition at line 5130 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::apply ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  X,
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
overridevirtual

Compute a sparse matrix-MultiVector multiply.

This method computes Y := beta*Y + alpha*Op(A)*X, where Op(A) is either $A$, $A^T$ (the transpose), or $A^H$ (the conjugate transpose).

If beta == 0, this operation will enjoy overwrite semantics: Y's entries will be ignored, and Y will be overwritten with the result of the multiplication, even if it contains NaN (not-a-number) floating-point entries.

Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 5091 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::hasTransposeApply ( ) const
overridevirtual

Whether apply() allows applying the transpose or conjugate transpose.

Reimplemented from Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 1079 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getDomainMap ( ) const
overridevirtual

The domain Map of this matrix.

This method implements Tpetra::Operator. If fillComplete() has not yet been called at least once on this matrix, or if the matrix was not constructed with a domain Map, then this method returns Teuchos::null.

Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 951 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRangeMap ( ) const
overridevirtual

The range Map of this matrix.

This method implements Tpetra::Operator. If fillComplete() has not yet been called at least once on this matrix, or if the matrix was not constructed with a domain Map, then this method returns Teuchos::null.

Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 958 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP< RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::add ( const Scalar &  alpha,
const RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Scalar &  beta,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  domainMap,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
) const
overridevirtual

Implementation of RowMatrix::add: return alpha*A + beta*this.

This override of the default implementation ensures that, when called on a CrsMatrix, this method always returns a CrsMatrix of exactly the same type as *this. "Exactly the same type" means that all the template parameters match, including the fifth template parameter. The input matrix A need not necessarily be a CrsMatrix or a CrsMatrix of the same type as *this, though this method may be able to optimize further in that case.

Reimplemented from Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 7401 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
std::string Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::description ( ) const
overridevirtual

A one-line description of this object.

Reimplemented from Tpetra::DistObject< char, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 5198 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
overridevirtual

Print this object with the given verbosity level to the given output stream.

Reimplemented from Tpetra::DistObject< char, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 5224 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::checkSizes ( const SrcDistObject source)
overridevirtual

Compare the source and target (this) objects for compatibility.

Returns
True if they are compatible, else false.

Implements Tpetra::DistObject< char, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 5460 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::copyAndPermute ( const SrcDistObject source,
const size_t  numSameIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteToLIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteFromLIDs,
const CombineMode  CM 
)
overrideprotectedvirtual

DistObject copyAndPermute has multiple overloads – use copyAndPermutes for anything we don't override

Definition at line 6018 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::unpackAndCombine ( const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  importLIDs,
Kokkos::DualView< char *, buffer_device_type imports,
Kokkos::DualView< size_t *, buffer_device_type numPacketsPerLID,
const size_t  constantNumPackets,
const CombineMode  CM 
)
override

Unpack the imported column indices and values, and combine into matrix.

Warning
The allowed CombineMode depends on whether the matrix's graph is static or dynamic. ADD, REPLACE, and ABSMAX are valid for a static graph, but INSERT is not. ADD and INSERT are valid for a dynamic graph; ABSMAX and REPLACE have not yet been implemented (and would require serious changes to matrix assembly in order to implement sensibly).

Definition at line 6871 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::packNew ( const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  exportLIDs,
Kokkos::DualView< char *, buffer_device_type > &  exports,
const Kokkos::DualView< size_t *, buffer_device_type > &  numPacketsPerLID,
size_t &  constantNumPackets 
) const

Pack this object's data for an Import or Export.

< DistObject has overloaded unpackAndCombine, use the DistObject's implementation for anything we don't override.

Warning
To be called only by the packAndPrepare method of appropriate classes of DistObject.
Parameters
exportLIDs[in] Local indices of the rows to pack.
exports[out] On output: array of packed matrix entries; allocated by method.
numPacketsPerLID[out] On output: numPacketsPerLID[i] is the number of bytes of the exports array used for storing packed local row exportLIDs[i].
constantNumPackets[out] If zero on output, the packed rows may have different numbers of entries. If nonzero on output, then that number gives the constant number of entries for all packed rows on all processes in the matrix's communicator.

The number of "packets" per row is the number of bytes per row. Each row has the following storage format:

[numEnt, vals, inds],

where:

  • numEnt (LocalOrdinal): number of entries in the row.
  • vals: array of Scalar. For the k-th entry in the row, vals[k] is its value and inds[k] its global column index.
  • inds: array of GlobalOrdinal. For the k-th entry in the row, vals[k] is its value and inds[k] its global column index.

We reserve the right to pad for alignment in the future. In that case, the number of bytes reported by numPacketsPerLID will reflect padding to align each datum to its size, and the row will have final padding as well to ensure that the next row is aligned. Rows with zero entries will still take zero bytes, however.

RowMatrix::pack will always use the same packing scheme as this method. This ensures correct Import / Export from a RowMatrix to a CrsMatrix.

We do not recommend relying on the details of this packing scheme. We describe it here more for Tpetra developers and less for users.

DistObject requires packing an object's entries as type Packet, which is the first template parameter of DistObject. Since sparse matrices have both values and indices, we use Packet=char and pack them into buffers of char (really "byte"). Indices are stored as global indices, in case the source and target matrices have different column Maps (or don't have a column Map yet).

Currently, we only pack values and column indices. Row indices are stored implicitly as the local indices (LIDs) to pack (see exportLIDs). This is because a DistObject instance only has one Map, and currently we use the row Map for CrsMatrix (and RowMatrix). This makes redistribution of matrices with 2-D distributions less efficient, but it works for now. This may change in the future.

On output, numPacketsPerLID[i] gives the number of bytes used to pack local row exportLIDs[i] of this object (the source object of an Import or Export). If offset is the exclusive prefix sum-scan of numPacketsPerLID, then on output, exports[offset[i] .. offset[i+1]] (half-exclusive range) contains the packed entries for local row exportLIDs[i].

Entries for each row use a "struct of arrays" pattern to match how sparse matrices actually store their data. The number of entries in the row goes first, all values go next, and all column indices (stored as global indices) go last. Values and column indices occur in the same order. Rows with zero entries always take zero bytes (we do not store their number of entries explicitly). This ensures sparsity of storage and communication in case most rows are empty.

GCC >= 4.9 and recent-future versions of the Intel compiler implement stricter aliasing rules that forbid unaligned type punning. If we were to pack as an "array of structs" – in this case, an array of (Scalar, GlobalOrdinal) pairs – then we would either have to pad each matrix entry for alignment, or call memcpy twice per matrix entry to pack and unpack. The "struct of arrays" storage scheme reduces the padding requirement to a constant per row, or reduces the number of memcpy calls to two per row.

We include the number of entries in each row in that row's packed data, to make unpacking easier. This saves us from an error-prone computation to find the number of entries from the number of bytes. That computation gets even more difficult if we have to introduce padding for alignment in the future. Knowing the number of entries for each row also makes parallelizing packing and unpacking easier.

Definition at line 6569 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
local_matrix_host_type::values_type::const_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalValuesHost ( Access::ReadOnlyStruct  s) const
inline

Get the Kokkos local values on host, read only.

Definition at line 3318 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
local_matrix_host_type::values_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalValuesHost ( Access::ReadWriteStruct  s)
inline

Get the Kokkos local values on host, read write.

Definition at line 3325 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
local_matrix_host_type::values_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalValuesHost ( Access::OverwriteAllStruct  s)
inline

Get the Kokkos local values on host, overwrite all.

Definition at line 3332 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
local_matrix_device_type::values_type::const_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalValuesDevice ( Access::ReadOnlyStruct  s) const
inline

Get the Kokkos local values on device, read only.

Definition at line 3339 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
local_matrix_device_type::values_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalValuesDevice ( Access::ReadWriteStruct  s)
inline

Get the Kokkos local values on device, read write.

Definition at line 3346 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
local_matrix_device_type::values_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalValuesDevice ( Access::OverwriteAllStruct  s)
inline

Get the Kokkos local values on device, overwrite all.

Definition at line 3353 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::importAndFillComplete ( Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &  destMatrix,
const import_type importer,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
) const

Import from this to the given destination matrix, and make the result fill complete.

If destMatrix.is_null(), this creates a new matrix as the destination. (This is why destMatrix is passed in by nonconst reference to RCP.) Otherwise it checks for "pristine" status and throws if that is not the case. "Pristine" means that the matrix has no entries and is not fill complete.

Use of the "non-member constructor" version of this method, exportAndFillCompleteCrsMatrix, is preferred for user applications.

Warning
This method is intended for expert developer use only, and should never be called by user code.

Definition at line 9121 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::importAndFillComplete ( Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &  destMatrix,
const import_type rowImporter,
const import_type domainImporter,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
) const

Import from this to the given destination matrix, and make the result fill complete.

If destMatrix.is_null(), this creates a new matrix as the destination. (This is why destMatrix is passed in by nonconst reference to RCP.) Otherwise it checks for "pristine" status and throws if that is not the case. "Pristine" means that the matrix has no entries and is not fill complete.

Use of the "non-member constructor" version of this method, exportAndFillCompleteCrsMatrix, is preferred for user applications.

Warning
This method is intended for expert developer use only, and should never be called by user code.

Definition at line 9133 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::exportAndFillComplete ( Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &  destMatrix,
const export_type exporter,
const Teuchos::RCP< const map_type > &  domainMap = Teuchos::null,
const Teuchos::RCP< const map_type > &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
) const

Export from this to the given destination matrix, and make the result fill complete.

If destMatrix.is_null(), this creates a new matrix as the destination. (This is why destMatrix is passed in by nonconst reference to RCP.) Otherwise it checks for "pristine" status and throws if that is not the case. "Pristine" means that the matrix has no entries and is not fill complete.

Use of the "non-member constructor" version of this method, exportAndFillCompleteCrsMatrix, is preferred for user applications.

Warning
This method is intended for expert developer use only, and should never be called by user code.

Definition at line 9146 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::exportAndFillComplete ( Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &  destMatrix,
const export_type rowExporter,
const export_type domainExporter,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
) const

Export from this to the given destination matrix, and make the result fill complete.

If destMatrix.is_null(), this creates a new matrix as the destination. (This is why destMatrix is passed in by nonconst reference to RCP.) Otherwise it checks for "pristine" status and throws if that is not the case. "Pristine" means that the matrix has no entries and is not fill complete.

Use of the "non-member constructor" version of this method, exportAndFillCompleteCrsMatrix, is preferred for user applications.

Warning
This method is intended for expert developer use only, and should never be called by user code.

Definition at line 9158 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal, class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::insertGlobalValuesImpl ( crs_graph_type graph,
RowInfo rowInfo,
const GlobalOrdinal  gblColInds[],
const impl_scalar_type  vals[],
const size_t  numInputEnt 
)
protectedvirtual

Common implementation detail of insertGlobalValues and insertGlobalValuesFiltered.

Precondition
&graph == this->getCrsGraph ().getRawPtr ()
rowInfo == graph.getRowInfo (rowInfo.localRow)
graph.getRowMap ()->isNodeLocalElement (rowInfo.localRow)
! this->isStaticGraph ()
If graph has a column Map, then all entries of gblColInds are in the column Map on the calling process. That is, the entries of gblColInds (and their corresponding vals entries) are "prefiltered," if we needed to filter them.

Definition at line 1886 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::allocateValues ( ELocalGlobal  lg,
GraphAllocationStatus  gas,
const bool  verbose 
)
protected

Allocate values (and optionally indices) using the Node.

Parameters
gas[in] If GraphNotYetAllocated, allocate the indices of myGraph_ via allocateIndices(lg) before allocating values.
lg[in] Argument passed into myGraph_->allocateIndices(), if applicable.
verbose[in] Whether to print verbose debugging output.
Precondition
If the graph (that is, staticGraph_) indices are already allocated, then gas must be GraphAlreadyAllocated. Otherwise, gas must be GraphNotYetAllocated. We only check for this precondition in debug mode.
If the graph indices are not already allocated, then the graph must be owned by the matrix.

Definition at line 1093 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mergeRowIndicesAndValues ( size_t  rowLen,
local_ordinal_type cols,
impl_scalar_type vals 
)
staticprotected

Merge duplicate row indices in the given row, along with their corresponding values.

This method is only called by sortAndMergeIndicesAndValues(), and only when the matrix owns the graph, not when the matrix was constructed with a const graph.

Precondition
The graph is not already storage optimized: isStorageOptimized() == false
Returns
The new row length, after merging.

Definition at line 4611 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::sortAndMergeIndicesAndValues ( const bool  sorted,
const bool  merged 
)
protected

Sort and merge duplicate local column indices in all rows on the calling process, along with their corresponding values.

Precondition
The matrix is locally indexed (more precisely, not globally indexed).
The matrix owns its graph.
The matrix's graph is not already storage optimized: isStorageOptimized() == false.
Parameters
sorted[in] If true, the column indices in each row on the calling process are already sorted.
merged[in] If true, the column indices in each row on the calling process are already merged.

Definition at line 4647 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::haveGlobalConstants ( ) const

Returns true if globalConstants have been computed; false otherwise.

Definition at line 4316 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getColumnMapMultiVector ( const MV X_domainMap,
const bool  force = false 
) const
protected

Create a (or fetch a cached) column Map MultiVector.

Parameters
X_domainMap[in] A domain Map Multivector. The returned MultiVector, if nonnull, will have the same number of columns as Y_domainMap.
force[in] Force creating the MultiVector if it hasn't been created already.

The force parameter is helpful when the domain Map and the column Map are the same (so that normally we wouldn't need the column Map MultiVector), but the following (for example) holds:

  1. The kernel needs a constant stride input MultiVector, but the given input MultiVector is not constant stride.

We don't test for the above in this method, because it depends on the specific kernel.

Definition at line 7275 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRowMapMultiVector ( const MV Y_rangeMap,
const bool  force = false 
) const
protected

Create a (or fetch a cached) row Map MultiVector.

Parameters
Y_rangeMap[in] A range Map Multivector. The returned MultiVector, if nonnull, will have the same number of columns as Y_rangeMap.
force[in] Force creating the MultiVector if it hasn't been created already.

The force parameter is helpful when the range Map and the row Map are the same (so that normally we wouldn't need the row Map MultiVector), but one of the following holds:

  1. The kernel needs a constant stride output MultiVector, but the given output MultiVector is not constant stride.
  2. The kernel does not permit aliasing of its input and output MultiVector arguments, but they do alias each other.

We don't test for the above in this method, because it depends on the specific kernel.

Definition at line 7332 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::applyNonTranspose ( const MV X_in,
MV Y_in,
Scalar  alpha,
Scalar  beta 
) const
protected

Special case of apply() for mode == Teuchos::NO_TRANS.

Definition at line 4707 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::applyTranspose ( const MV X_in,
MV Y_in,
const Teuchos::ETransp  mode,
Scalar  alpha,
Scalar  beta 
) const
protected

Special case of apply() for mode != Teuchos::NO_TRANS.

Definition at line 4874 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::values_dualv_type::t_host::const_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getValuesViewHost ( const RowInfo rowinfo) const
protected

Get a const Host view of the locally owned values row myRow, such that rowinfo = getRowInfo(myRow).

Definition at line 3133 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::values_dualv_type::t_dev::const_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getValuesViewDevice ( const RowInfo rowinfo) const
protected

Get a const Device view of the locally owned values row myRow, such that rowinfo = getRowInfo(myRow).

Definition at line 3161 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::values_dualv_type::t_host Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getValuesViewHostNonConst ( const RowInfo rowinfo)
protected

Get a non-const Host view of the locally owned values row myRow, such that rowinfo = getRowInfo(myRow).

Definition at line 3147 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::values_dualv_type::t_dev Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getValuesViewDeviceNonConst ( const RowInfo rowinfo)
protected

Get a non-const Device view of the locally owned values row myRow, such that rowinfo = getRowInfo(myRow).

Definition at line 3175 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::swap ( CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  matrix)
protected

Swaps the data from *this with the data and maps from crsMatrix.

Parameters
matrix[in/out] a crsMatrix

Definition at line 793 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fillLocalMatrix ( const Teuchos::RCP< Teuchos::ParameterList > &  params)
protected

Fill data into the local matrix.

This method is only called in fillComplete(), and it is only called if the graph's structure is already fixed (that is, if the matrix does not own the graph).

Definition at line 1582 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fillLocalGraphAndMatrix ( const Teuchos::RCP< Teuchos::ParameterList > &  params)
protected

Fill data into the local graph and matrix.

This method is only called in fillComplete(), and it is only called if the graph's structure is not already fixed (that is, if the matrix does own the graph).

Definition at line 1203 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::checkInternalState ( ) const
protected

Check that this object's state is sane; throw if it's not.

Definition at line 5161 of file Tpetra_CrsMatrix_def.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowCopy ( GlobalOrdinal  GlobalRow,
nonconst_global_inds_host_view_type &  Indices,
nonconst_values_host_view_type &  Values,
size_t &  NumEntries 
) const
pure virtualinherited

Get a copy of the given global row's entries.

This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Parameters
GlobalRow[in] Global index of the row.
Indices[out] Global indices of the columns corresponding to values.
Values[out] Matrix values.
NumEntries[out] Number of stored entries on the calling process; length of Indices and Values.

This method throws std::runtime_error if either Indices or Values is not large enough to hold the data associated with row GlobalRow. If GlobalRow does not belong to the calling process, then the method sets NumIndices to Teuchos::OrdinalTraits<size_t>::invalid(), and does not modify Indices or Values.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowCopy ( LocalOrdinal  LocalRow,
nonconst_local_inds_host_view_type &  Indices,
nonconst_values_host_view_type &  Values,
size_t &  NumEntries 
) const
pure virtualinherited

Get a copy of the given local row's entries.

This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Parameters
LocalRow[in] Local index of the row.
Indices[out] Local indices of the columns corresponding to values.
Values[out] Matrix values.
NumEntries[out] Number of stored entries on the calling process; length of Indices and Values.

This method throws std::runtime_error if either Indices or Values is not large enough to hold the data associated with row LocalRow. If LocalRow does not belong to the calling process, then the method sets NumIndices to Teuchos::OrdinalTraits<size_t>::invalid(), and does not modify Indices or Values.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowView ( GlobalOrdinal  GlobalRow,
global_inds_host_view_type &  indices,
values_host_view_type &  values 
) const
pure virtualinherited

Get a constant, nonpersisting, globally indexed view of the given row of the matrix.

The returned views of the column indices and values are not guaranteed to persist beyond the lifetime of this. Furthermore, some RowMatrix implementations allow changing the values, or the indices and values. Any such changes invalidate the returned views.

This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Precondition
isGloballyIndexed () && supportsRowViews ()
Postcondition
indices.size () == getNumEntriesInGlobalRow (GlobalRow)
Parameters
GlobalRow[in] Global index of the row.
Indices[out] Global indices of the columns corresponding to values.
Values[out] Matrix values.

If GlobalRow does not belong to this node, then indices is set to null.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getLocalRowView ( LocalOrdinal  LocalRow,
local_inds_host_view_type &  indices,
values_host_view_type &  values 
) const
pure virtualinherited

Get a constant, nonpersisting, locally indexed view of the given row of the matrix.

The returned views of the column indices and values are not guaranteed to persist beyond the lifetime of this. Furthermore, some RowMatrix implementations allow changing the values, or the indices and values. Any such changes invalidate the returned views.

This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Precondition
isLocallyIndexed () && supportsRowViews ()
Postcondition
indices.size () == getNumEntriesInGlobalRow (LocalRow)
Parameters
LocalRow[in] Local index of the row.
Indices[out] Local indices of the columns corresponding to values.
Values[out] Matrix values.

If LocalRow does not belong to this node, then indices is set to null.

template<class Scalar , class LocalOrdinal, class GlobalOrdinal , class Node >
void Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::pack ( const Teuchos::ArrayView< const LocalOrdinal > &  exportLIDs,
Teuchos::Array< char > &  exports,
const Teuchos::ArrayView< size_t > &  numPacketsPerLID,
size_t &  constantNumPackets 
) const
virtualinherited

Pack this object's data for an Import or Export.

Warning
To be called only by the packAndPrepare method of appropriate classes of DistObject.

Subclasses may override this method to speed up or otherwise improve the implementation by exploiting more specific details of the subclass.

Implements Tpetra::Packable< char, LocalOrdinal >.

Definition at line 301 of file Tpetra_RowMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >::hasDiagonal ( ) const
virtualinherited

Whether this operator can return its diagonal.

By default, this returns false. Subclasses must override this method if they can supply a diagonal.

Definition at line 156 of file Tpetra_Operator.hpp.

void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::doImport ( const SrcDistObject< char, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Import< LocalOrdinal, GlobalOrdinal, Node > &  importer,
const CombineMode  CM,
const bool  restrictedMode = false 
)
inherited

Import data into this object using an Import object ("forward mode").

The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.

If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use this method with your precomputed Import object if you want to do an Import, else use doExport() with a precomputed Export object.

"Restricted Mode" does two things:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *this, in a "locallyFitted" sense.

This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.

Parameters
source[in] The "source" object for redistribution.
importer[in] Precomputed data redistribution plan. Its source Map must be the same as the input DistObject's Map, and its target Map must be the same as this->getMap().
CM[in] How to combine incoming data with the same global index.
void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::doImport ( const SrcDistObject< char, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Export< LocalOrdinal, GlobalOrdinal, Node > &  exporter,
const CombineMode  CM,
const bool  restrictedMode = false 
)
inherited

Import data into this object using an Export object ("reverse mode").

The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.

If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use the version of doImport() that takes a precomputed Import object in that case.

"Restricted Mode" does two things:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *this, in a "locallyFitted" sense.

This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.

Parameters
source[in] The "source" object for redistribution.
exporter[in] Precomputed data redistribution plan. Its target Map must be the same as the input DistObject's Map, and its source Map must be the same as this->getMap(). (Note the difference from forward mode.)
CM[in] How to combine incoming data with the same global index.
void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::doExport ( const SrcDistObject< char, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Export< LocalOrdinal, GlobalOrdinal, Node > &  exporter,
const CombineMode  CM,
const bool  restrictedMode = false 
)
inherited

Export data into this object using an Export object ("forward mode").

The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.

If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use this method with your precomputed Export object if you want to do an Export, else use doImport() with a precomputed Import object.

"Restricted Mode" does two things:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *this, in a "locallyFitted" sense.

This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.

Parameters
source[in] The "source" object for redistribution.
exporter[in] Precomputed data redistribution plan. Its source Map must be the same as the input DistObject's Map, and its target Map must be the same as this->getMap().
CM[in] How to combine incoming data with the same global index.
void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::doExport ( const SrcDistObject< char, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Import< LocalOrdinal, GlobalOrdinal, Node > &  importer,
const CombineMode  CM,
const bool  restrictedMode = false 
)
inherited

Export data into this object using an Import object ("reverse mode").

The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.

If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use the version of doExport() that takes a precomputed Export object in that case.

"Restricted Mode" does two things:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *this, in a "locallyFitted" sense.

This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.

Parameters
source[in] The "source" object for redistribution.
importer[in] Precomputed data redistribution plan. Its target Map must be the same as the input DistObject's Map, and its source Map must be the same as this->getMap(). (Note the difference from forward mode.)
CM[in] How to combine incoming data with the same global index.
bool Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::transferArrived ( ) const
inherited

Whether the data from an import/export operation has arrived, and is ready for the unpack and combine step.

bool Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::isDistributed ( ) const
inherited

Whether this is a globally distributed object.

For a definition of "globally distributed" (and its opposite, "locally replicated"), see the documentation of Map's isDistributed() method.

virtual Teuchos::RCP<const map_type> Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::getMap ( ) const
inlinevirtualinherited

The Map describing the parallel distribution of this object.

Note that some Tpetra objects might be distributed using multiple Map objects. For example, CrsMatrix has both a row Map and a column Map. It is up to the subclass to decide which Map to use when invoking the DistObject constructor.

Definition at line 589 of file Tpetra_DistObject_decl.hpp.

void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::print ( std::ostream &  os) const
inherited

Print this object to the given output stream.

We generally assume that all MPI processes can print to the given stream.

virtual void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::removeEmptyProcessesInPlace ( const Teuchos::RCP< const map_type > &  newMap)
virtualinherited

Remove processes which contain no entries in this object's Map.

Warning
This method is ONLY for use by experts. We highly recommend using the nonmember function of the same name defined in this file.
We make NO promises of backwards compatibility. This method may change or disappear at any time.

On input, this object is distributed over the Map returned by getMap() (the "original Map," with its communicator, the "original communicator"). The input newMap of this method must be the same as the result of calling getMap()->removeEmptyProcesses(). On processes in the original communicator which contain zero entries ("excluded processes," as opposed to "included processes"), the input newMap must be Teuchos::null (which is what getMap()->removeEmptyProcesses() returns anyway).

On included processes, reassign this object's Map (that would be returned by getMap()) to the input newMap, and do any work that needs to be done to restore correct semantics. On excluded processes, free any data that needs freeing, and do any other work that needs to be done to restore correct semantics.

This method has collective semantics over the original communicator. On exit, the only method of this object which is safe to call on excluded processes is the destructor. This implies that subclasses' destructors must not contain communication operations.

Returns
The object's new Map. Its communicator is a new communicator, distinct from the old Map's communicator, which contains a subset of the processes in the old communicator.
Note
The name differs from Map's method removeEmptyProcesses(), in order to emphasize that the operation on DistObject happens in place, modifying the input, whereas the operation removeEmptyProcess() on Map does not modify the input.
To implementers of DistObject subclasses: The default implementation of this class throws std::logic_error.
virtual size_t Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::constantNumberOfPackets ( ) const
protectedvirtualinherited

Whether the implementation's instance promises always to have a constant number of packets per LID (local index), and if so, how many packets per LID there are.

If this method returns zero, the instance says that it might possibly have a different number of packets for each LID (local index) to send or receive. If it returns nonzero, the instance promises that the number of packets is the same for all LIDs, and that the return value is this number of packets per LID.

The default implementation of this method returns zero. This does not affect the behavior of doTransfer() in any way. If a nondefault implementation returns nonzero, doTransfer() will use this information to avoid unnecessary allocation and / or resizing of arrays.

virtual void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::doTransfer ( const SrcDistObject< char, LocalOrdinal, GlobalOrdinal, Node > &  src,
const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &  transfer,
const char  modeString[],
const ReverseOption  revOp,
const CombineMode  CM,
const bool  restrictedMode 
)
protectedvirtualinherited

Redistribute data across (MPI) processes.

Parameters
src[in] The source object, to redistribute into the target object, which is *this object.
transfer[in] The Export or Import object representing the communication pattern. (Details::Transfer is the common base class of these two objects.)
modeString[in] Human-readable string, for verbose debugging output and error output, explaining what function called this method. Example: "doImport (forward)", "doExport (reverse)".
revOp[in] Whether to do a forward or reverse mode redistribution.
CM[in] The combine mode that describes how to combine values that map to the same global ID on the same process.
virtual bool Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::reallocArraysForNumPacketsPerLid ( const size_t  numExportLIDs,
const size_t  numImportLIDs 
)
protectedvirtualinherited

Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary.

Parameters
numExportLIDs[in] Number of entries in the exportLIDs input array argument of doTransfer().
numImportLIDs[in] Number of entries in the remoteLIDs input array argument of doTransfer().
Returns
Whether we actually reallocated either of the arrays.
Warning
This is an implementation detail of doTransferNew(). This needs to be protected, but that doesn't mean users should call this method.
void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::beginTransfer ( const SrcDistObject< char, LocalOrdinal, GlobalOrdinal, Node > &  src,
const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &  transfer,
const char  modeString[],
const ReverseOption  revOp,
const CombineMode  CM,
const bool  restrictedMode 
)
protectedinherited

Implementation detail of doTransfer.

LID DualViews come from the Transfer object given to doTransfer. They are always sync'd on both host and device. Users must never attempt to modify or sync them.

virtual void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::copyAndPermute ( const SrcDistObject< char, LocalOrdinal, GlobalOrdinal, Node > &  source,
const size_t  numSameIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteToLIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteFromLIDs,
const CombineMode  CM 
)
protectedvirtualinherited

Perform copies and permutations that are local to the calling (MPI) process.

Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this, packs the source object's data.

Precondition
permuteToLIDs and permuteFromLIDs are sync'd to both host and device. That is, permuteToLIDs.need_sync_host(), permuteToLIDs.need_sync_device(), permuteFromLIDs.need_sync_host(), and permuteFromLIDs.need_sync_device() are all false.
Parameters
source[in] On entry, the source object of the Export or Import operation.
numSameIDs[in] The number of elements that are the same on the source and target objects. These elements live on the same process in both the source and target objects.
permuteToLIDs[in] List of the elements that are permuted. They are listed by their local index (LID) in the destination object.
permuteFromLIDs[in] List of the elements that are permuted. They are listed by their local index (LID) in the source object.
CM[in] CombineMode to be used during copyAndPermute; may or may not be used by the particular object being called; behavior with respect to CombineMode may differ by object.
virtual void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::copyAndPermute ( const SrcDistObject< char, LocalOrdinal, GlobalOrdinal, Node > &  source,
const size_t  numSameIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteToLIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteFromLIDs,
const CombineMode  CM,
const execution_space space 
)
protectedvirtualinherited

Same as copyAndPermute, but do operations in space.

virtual void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::packAndPrepare ( const SrcDistObject< char, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  exportLIDs,
Kokkos::DualView< packet_type *, buffer_device_type > &  exports,
Kokkos::DualView< size_t *, buffer_device_type numPacketsPerLID,
size_t &  constantNumPackets 
)
protectedvirtualinherited

Pack data and metadata for communication (sends).

Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this, packs the source object's data.

Precondition
exportLIDs is sync'd to both host and device. That is, exportLIDs.need_sync_host () and exportLIDs.need_sync_device() are both false.
Parameters
source[in] Source object for the redistribution.
exportLIDs[in] List of the entries (as local IDs in the source object) that Tpetra will send to other processes.
exports[out] On exit, the packed data to send. Implementations must reallocate this as needed (prefer reusing the existing allocation if possible), and may modify and/or sync this wherever they like.
numPacketsPerLID[out] On exit, the implementation of this method must do one of two things: either set numPacketsPerLID[i] to the number of packets to be packed for exportLIDs[i] and set constantNumPackets to zero, or set constantNumPackets to a nonzero value. If the latter, the implementation must not modify the entries of numPacketsPerLID. If the former, the implementation may sync numPacketsPerLID this wherever it likes, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference.
constantNumPackets[out] On exit, 0 if the number of packets per LID could differ, else (if nonzero) the number of packets per LID (which must be constant).
virtual void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::packAndPrepare ( const SrcDistObject< char, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  exportLIDs,
Kokkos::DualView< packet_type *, buffer_device_type > &  exports,
Kokkos::DualView< size_t *, buffer_device_type numPacketsPerLID,
size_t &  constantNumPackets,
const execution_space space 
)
protectedvirtualinherited

Same as packAndPrepare, but in an execution space instance.

virtual void Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::unpackAndCombine ( const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  importLIDs,
Kokkos::DualView< packet_type *, buffer_device_type imports,
Kokkos::DualView< size_t *, buffer_device_type numPacketsPerLID,
const size_t  constantNumPackets,
const CombineMode  combineMode 
)
protectedvirtualinherited

Perform any unpacking and combining after communication.

Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this, unpacks the received data into itself, possibly modifying its entries.

Precondition
importLIDs is sync'd to both host and device. That is, importLIDs.need_sync_host () and importLIDs.need_sync_device() are both false.
Parameters
importLIDs[in] List of the entries (as LIDs in the destination object) we received from other processes.
imports[in/out] On input: Buffer of received data to unpack. DistObject promises nothing about where this is sync'd. Implementations may sync this wherever they like, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference.
numPacketsPerLID[in/out] On input: If constantNumPackets is zero, then numPacketsPerLID[i] contains the number of packets imported for importLIDs[i]. DistObject promises nothing about where this is sync'd. Implementations may sync this wherever they like, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference.
constantNumPackets[in] If nonzero, then the number of packets per LID is the same for all entries ("constant") and constantNumPackets is that number. If zero, then numPacketsPerLID[i] is the number of packets to unpack for LID importLIDs[i].
combineMode[in] The CombineMode to use when combining the imported entries with existing entries.
virtual bool Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::reallocImportsIfNeeded ( const size_t  newSize,
const bool  verbose,
const std::string *  prefix,
const bool  remoteLIDsContiguous = false,
const CombineMode  CM = INSERT 
)
protectedvirtualinherited

Reallocate imports_ if needed.

This unfortunately must be declared protected, for the same reason that imports_ is declared protected.

Parameters
newSize[in] New size of imports_.
verbose[in] Whether to print verbose debugging output to stderr on every (MPI) process in the communicator.
prefix[in] If verbose is true, then this is a nonnull prefix to print at the beginning of each line of verbose debugging output. Otherwise, not used.
Returns
Whether we actually reallocated.

We don't need a "reallocExportsIfNeeded" method, because exports_ always gets passed into packAndPrepare() by nonconst reference. Thus, that method can resize the DualView without needing to call other DistObject methods.

Friends And Related Function Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
template<class MatrixArray , class MultiVectorArray >
void batchedApply ( const MatrixArray &  Matrices,
const typename std::remove_pointer< typename MultiVectorArray::value_type >::type &  X,
MultiVectorArray &  Y,
typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type  alpha = Teuchos::ScalarTraits< typename std::remove_pointer<typename MatrixArray::value_type>::type::scalar_type>::one(),
typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type  beta = Teuchos::ScalarTraits< typename std::remove_pointer<typename MatrixArray::value_type>::type::scalar_type>::zero(),
Teuchos::RCP< Teuchos::ParameterList >  params = Teuchos::null 
)
friend

Does multiply matrix apply() calls with a single X vector.

Computes Y[i] = beta*Y[i] + alpha * Matrices[i] * X, for i=1,...,N

This routine only communicates the interprocessor portions of X once, if such a thing is possible (aka the ColMap's of the matrices match).

If the Y's are replicated this will not use the reduced communication path.

If X is aliased to any Y but the last one, this routine will throw

Parameters
Matrices[in] - [std::vector|Teuchos::Array|Teuchos::ArrayRCP] of Tpetra::CrsMatrix objects. These matrices can different numbers of rows.
X[in] - Tpetra::MultiVector or Tpetra::Vector object.
Y[out] - [std::vector|Teuchos::Array|Teuchos::ArrayRCP] of Tpetra::MultiVector or Tpetra::Vector objects. These must have the same number of vectors as X.
alpha[in] - alpha parameter. Defaults to one.
beta[in] - beta parameter. Defaults to zero.
params[in/out] - The "can batch" parameter can either be unset, be true or be false on input. If it is unset, maps will be checked with isSameAs() for compatibility. If true, the map check will be skipped and batching will be used if none of the cheap checks fail. If false, batching will not be used. This parameter will be set on output to either true or false depending on if batching was used during this call. Defaults to NULL.

Definition at line 74 of file Tpetra_Apply_Helpers.hpp.

template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > importAndFillCompleteCrsMatrix ( const Teuchos::RCP< const CrsMatrixType > &  sourceMatrix,
const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  importer,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  domainMap = Teuchos::null,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)
related

Nonmember CrsMatrix constructor that fuses Import and fillComplete().

Template Parameters
CrsMatrixTypeA specialization of CrsMatrix.

A common use case is to create an empty destination CrsMatrix, redistribute from a source CrsMatrix (by an Import or Export operation), then call fillComplete() on the destination CrsMatrix. This constructor fuses these three cases, for an Import redistribution.

Fusing redistribution and fillComplete() exposes potential optimizations. For example, it may make constructing the column Map faster, and it may avoid intermediate unoptimized storage in the destination CrsMatrix. These optimizations may improve performance for specialized kernels like sparse matrix-matrix multiply, as well as for redistributing data after doing load balancing.

The resulting matrix is fill complete (in the sense of isFillComplete()) and has optimized storage (in the sense of isStorageOptimized()). By default, its domain Map is the domain Map of the source matrix, and its range Map is the range Map of the source matrix.

Warning
If the target Map of the Import is a subset of the source Map of the Import, then you cannot use the default range Map. You should instead construct a nonoverlapping version of the target Map and supply that as the nondefault value of the range Map.
Parameters
sourceMatrix[in] The source matrix from which to import. The source of an Import must have a nonoverlapping distribution.
importer[in] The Import instance containing a precomputed redistribution plan. The source Map of the Import must be the same as the rowMap of sourceMatrix unless the "Reverse Mode" option on the params list, in which case the targetMap of Import must match the rowMap of the sourceMatrix
domainMap[in] Domain Map of the returned matrix. If null, we use the default, which is the domain Map of the source matrix.
rangeMap[in] Range Map of the returned matrix. If null, we use the default, which is the range Map of the source matrix.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.
template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > importAndFillCompleteCrsMatrix ( const Teuchos::RCP< const CrsMatrixType > &  sourceMatrix,
const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  rowImporter,
const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  domainImporter,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  domainMap,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
)
related

Nonmember CrsMatrix constructor that fuses Import and fillComplete().

Template Parameters
CrsMatrixTypeA specialization of CrsMatrix.

A common use case is to create an empty destination CrsMatrix, redistribute from a source CrsMatrix (by an Import or Export operation), then call fillComplete() on the destination CrsMatrix. This constructor fuses these three cases, for an Import redistribution.

Fusing redistribution and fillComplete() exposes potential optimizations. For example, it may make constructing the column Map faster, and it may avoid intermediate unoptimized storage in the destination CrsMatrix. These optimizations may improve performance for specialized kernels like sparse matrix-matrix multiply, as well as for redistributing data after doing load balancing.

The resulting matrix is fill complete (in the sense of isFillComplete()) and has optimized storage (in the sense of isStorageOptimized()). By default, its domain Map is the domain Map of the source matrix, and its range Map is the range Map of the source matrix.

Warning
If the target Map of the Import is a subset of the source Map of the Import, then you cannot use the default range Map. You should instead construct a nonoverlapping version of the target Map and supply that as the nondefault value of the range Map.
Parameters
sourceMatrix[in] The source matrix from which to import. The source of an Import must have a nonoverlapping distribution.
rowImporter[in] The Import instance containing a precomputed redistribution plan. The source Map of the Import must be the same as the rowMap of sourceMatrix unless the "Reverse Mode" option on the params list, in which case the targetMap of Import must match the rowMap of the sourceMatrix
domainImporter[in] The Import instance containing a precomputed redistribution plan. The source Map of the Import must be the same as the domainMap of sourceMatrix unless the "Reverse Mode" option on the params list, in which case the targetMap of Import must match the domainMap of the sourceMatrix
domainMap[in] Domain Map of the returned matrix.
rangeMap[in] Range Map of the returned matrix.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.
template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > exportAndFillCompleteCrsMatrix ( const Teuchos::RCP< const CrsMatrixType > &  sourceMatrix,
const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  exporter,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  domainMap = Teuchos::null,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)
related

Nonmember CrsMatrix constructor that fuses Export and fillComplete().

Template Parameters
CrsMatrixTypeA specialization of CrsMatrix.

For justification, see the documentation of importAndFillCompleteCrsMatrix() (which is the Import analog of this function).

The resulting matrix is fill complete (in the sense of isFillComplete()) and has optimized storage (in the sense of isStorageOptimized()). By default, its domain Map is the domain Map of the source matrix, and its range Map is the range Map of the source matrix.

Parameters
sourceMatrix[in] The source matrix from which to export. Its row Map may be overlapping, since the source of an Export may be overlapping.
exporter[in] The Export instance containing a precomputed redistribution plan. The source Map of the Export must be the same as the row Map of sourceMatrix.
domainMap[in] Domain Map of the returned matrix. If null, we use the default, which is the domain Map of the source matrix.
rangeMap[in] Range Map of the returned matrix. If null, we use the default, which is the range Map of the source matrix.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.
template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > exportAndFillCompleteCrsMatrix ( const Teuchos::RCP< const CrsMatrixType > &  sourceMatrix,
const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  rowExporter,
const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  domainExporter,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  domainMap,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
)
related

Nonmember CrsMatrix constructor that fuses Export and fillComplete().

Template Parameters
CrsMatrixTypeA specialization of CrsMatrix.

For justification, see the documentation of importAndFillCompleteCrsMatrix() (which is the Import analog of this function).

The resulting matrix is fill complete (in the sense of isFillComplete()) and has optimized storage (in the sense of isStorageOptimized()). By default, its domain Map is the domain Map of the source matrix, and its range Map is the range Map of the source matrix.

Parameters
sourceMatrix[in] The source matrix from which to export. Its row Map may be overlapping, since the source of an Export may be overlapping.
rowExporter[in] The Export instance containing a precomputed redistribution plan. The source Map of the Export must be the same as the row Map of sourceMatrix.
domainExporter[in] The Export instance containing a precomputed redistribution plan. The source Map of the Export must be the same as the domain Map of sourceMatrix.
domainMap[in] Domain Map of the returned matrix.
rangeMap[in] Range Map of the returned matrix.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > createCrsMatrix ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  map,
const size_t  maxNumEntriesPerRow = 0,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)
related

Create an empty CrsMatrix given a row map and a single integer upper bound on the number of stored entries per row.

Definition at line 4070 of file Tpetra_CrsMatrix_decl.hpp.

template<class CrsMatrixType >
void removeCrsMatrixZeros ( CrsMatrixType &  matrix,
typename Teuchos::ScalarTraits< typename CrsMatrixType::scalar_type >::magnitudeType const &  threshold = Teuchos::ScalarTraits<typename CrsMatrixType::scalar_type>::magnitude( Teuchos::ScalarTraits<typename CrsMatrixType::scalar_type>::zero() ) 
)
related

Remove zero entries from a matrix.

Parameters
matrix[in/out] CrsMatrix
threshold[in] magnitude threshold below which an entry is deemed to be zero

Definition at line 4171 of file Tpetra_CrsMatrix_decl.hpp.

Member Data Documentation

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
ordinal_rowptrs_type Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ordinalRowptrs
mutableprotected

local_ordinal typed version of local matrix's rowptrs. This allows the LocalCrsMatrixOperator to have rowptrs and entries be the same type, so cuSPARSE SpMV (including merge-path) can be used for apply. This is allocated and populated lazily in getLocalMultiplyOperator(), only if all 4 conditions are met:

  • node_type is KokkosCudaWrapperNode
  • the cuSPARSE TPL is enabled
  • local_ordinal_type can represent getLocalNumEntries()

Definition at line 2557 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP<MV> Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::importMV_
mutableprotected

Column Map MultiVector used in apply().

This is a column Map MultiVector. It is used as the target of the forward mode Import operation (if necessary) in apply(), and the source of the reverse mode Export operation (if necessary) in these methods. Both of these methods create this MultiVector on demand if needed, and reuse it (if possible) for subsequent calls.

This is declared mutable because the methods in question are const, yet want to cache the MultiVector for later use.

Definition at line 3830 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP<MV> Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::exportMV_
mutableprotected

Row Map MultiVector used in apply().

This is a row Map MultiVector. It is uses as the source of the forward mode Export operation (if necessary) in apply(), and the target of the reverse mode Import operation (if necessary) in these methods. Both of these methods create this MultiVector on demand if needed, and reuse it (if possible) for subsequent calls.

This is declared mutable because the methods in question are const, yet want to cache the MultiVector for later use.

Definition at line 3844 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
Details::EStorageStatus Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::storageStatus_
protected
Initial value:
=
Details::STORAGE_1D_UNPACKED

Status of the matrix's storage, when not in a fill-complete state.

The phrase "When not in a fill-complete state" is important. When the matrix is fill complete, it always uses 1-D "packed" storage. However, if the "Optimize Storage" parameter to fillComplete was false, the matrix may keep unpacked 1-D storage around and resume it on the next resumeFill call.

Definition at line 3988 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
bool Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fillComplete_ = false
protected

Whether the matrix is fill complete.

Definition at line 3992 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
std::map<GlobalOrdinal, std::pair<Teuchos::Array<GlobalOrdinal>, Teuchos::Array<Scalar> > > Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::nonlocals_
protected

Nonlocal data added using insertGlobalValues().

These data are cleared by globalAssemble(), once it finishes redistributing them to their owning processes.

For a given nonowned global row gRow which was given to insertGlobalValues() or sumIntoGlobalValues(), nonlocals_[gRow].first[k] is the column index of an inserted entry, and nonlocals_[gRow].second[k] is its value. Duplicate column indices for the same row index are allowed and will be summed during globalAssemble().

This used to be a map from GlobalOrdinal to (GlobalOrdinal, Scalar) pairs. This makes gcc issue a "note" about the ABI of structs containing std::complex members changing. CDash reports this as a warning, even though it's a "note," not a warning. However, I don't want it to show up, so I rearranged the map's value type to a pair of arrays, rather than an array of pairs.

Note
For Epetra developers: Tpetra::CrsMatrix corresponds more to Epetra_FECrsMatrix than to Epetra_CrsMatrix. The insertGlobalValues() method in Tpetra::CrsMatrix, unlike its corresponding method in Epetra_CrsMatrix, allows insertion into rows which are not owned by the calling process. The globalAssemble() method redistributes these to their owning processes.

Definition at line 4022 of file Tpetra_CrsMatrix_decl.hpp.

Teuchos::RCP<const map_type> Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::map_
protectedinherited

The Map over which this object is distributed.

Definition at line 998 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<packet_type*, buffer_device_type> Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::imports_
protectedinherited

Buffer into which packed data are imported (received from other processes).

Unfortunately, I had to declare these protected, because CrsMatrix uses them at one point. Please, nobody else use them.

Definition at line 1011 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<size_t*, buffer_device_type> Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::numImportPacketsPerLID_
protectedinherited

Number of packets to receive for each receive operation.

This array is used in Distributor::doPosts() (and doReversePosts()) when starting the ireceive operation.

This may be ignored in doTransfer() if constantNumPackets is nonzero, indicating a constant number of packets per LID. (For example, MultiVector sets the constantNumPackets output argument of packAndPrepare() to the number of columns in the multivector.)

Unfortunately, I had to declare this protected, because CrsMatrix uses it at one point. Please, nobody else use it.

Definition at line 1051 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<packet_type*, buffer_device_type> Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::exports_
protectedinherited

Buffer from which packed data are exported (sent to other processes).

Unfortunately, I had to declare this protected, because CrsMatrix uses it at one point. Please, nobody else use it.

Definition at line 1058 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<size_t*, buffer_device_type> Tpetra::DistObject< char , LocalOrdinal, GlobalOrdinal, Node >::numExportPacketsPerLID_
protectedinherited

Number of packets to send for each send operation.

This array is used in Distributor::doPosts() (and doReversePosts()) for preparing for the send operation.

This may be ignored in doTransfer() if constantNumPackets is nonzero, indicating a constant number of packets per LID. (For example, MultiVector sets the constantNumPackets output argument of packAndPrepare() to the number of columns in the multivector.)

Unfortunately, I had to declare this protected, because CrsMatrix uses them at one point. Please, nobody else use it.

Definition at line 1073 of file Tpetra_DistObject_decl.hpp.


The documentation for this class was generated from the following files: