|
Tpetra parallel linear algebra
Version of the Day
|
Go to the documentation of this file.
42 #ifndef TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_DECL_HPP
43 #define TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_DECL_HPP
48 #include "Tpetra_CrsGraph.hpp"
49 #include "Tpetra_RowMatrix.hpp"
50 #include "Tpetra_Experimental_BlockMultiVector.hpp"
54 namespace Experimental {
134 template<
class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type,
136 class GO = ::Tpetra::Details::DefaultTypes::global_ordinal_type,
145 using STS = Teuchos::ScalarTraits<Scalar>;
186 typedef ::Tpetra::MultiVector<Scalar, LO, GO, node_type>
mv_type;
194 Kokkos::MemoryTraits<Kokkos::Unmanaged> >
200 Kokkos::MemoryTraits<Kokkos::Unmanaged> >
206 Kokkos::MemoryTraits<Kokkos::Unmanaged> >
212 Kokkos::MemoryTraits<Kokkos::Unmanaged> >
259 Teuchos::RCP<const map_type>
getRowMap ()
const;
262 Teuchos::RCP<const map_type>
getColMap ()
const;
284 Teuchos::ETransp mode = Teuchos::NO_TRANS,
285 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one (),
286 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero ())
const;
330 describe (Teuchos::FancyOStream& out,
331 const Teuchos::EVerbosityLevel verbLevel)
const;
341 virtual Teuchos::RCP<const ::Tpetra::RowGraph<LO,GO,Node> >
getGraph ()
const;
352 Teuchos::ETransp mode = Teuchos::NO_TRANS,
353 const Scalar alpha = Teuchos::ScalarTraits<Scalar>::one (),
354 const Scalar beta = Teuchos::ScalarTraits<Scalar>::zero ());
361 const ::Tpetra::MultiVector<Scalar,LO,GO,Node> &B,
362 const ::Tpetra::MultiVector<Scalar,LO,GO,Node> &D,
363 const Scalar& dampingFactor,
366 const bool zeroInitialGuess)
const;
373 const ::Tpetra::MultiVector<Scalar,LO,GO,Node>& B,
374 const ::Tpetra::MultiVector<Scalar,LO,GO,Node>& D,
375 const Teuchos::ArrayView<LO>& rowIndices,
376 const Scalar& dampingFactor,
379 const bool zeroInitialGuess)
const;
404 Kokkos::MemoryUnmanaged>& D_inv,
438 const LO numColInds)
const;
470 const LO numColInds)
const;
510 Teuchos::ArrayView<const LO> &indices,
511 Teuchos::ArrayView<const Scalar> &values)
const;
516 const Teuchos::ArrayView<LO> &Indices,
517 const Teuchos::ArrayView<Scalar> &Values,
518 size_t &NumEntries)
const;
521 getLocalBlock (
const LO localRowInd,
const LO localColInd)
const;
550 const LO numColInds)
const;
559 const ptrdiff_t offsets[],
561 const LO numOffsets)
const;
570 const ptrdiff_t offsets[],
572 const LO numOffsets)
const;
617 return (*errs_).is_null () ? std::string (
"") : (*errs_)->str ();
653 Kokkos::MemoryUnmanaged>& offsets)
const;
660 void TPETRA_DEPRECATED
678 Kokkos::MemoryUnmanaged>& diag,
680 Kokkos::MemoryUnmanaged>& offsets)
const;
697 Kokkos::MemoryUnmanaged>& diag,
698 const Teuchos::ArrayView<const size_t>& offsets)
const;
713 void TPETRA_DEPRECATED
715 const Teuchos::ArrayView<const size_t>& offsets)
const;
723 const LO numColInds)
const;
728 const ptrdiff_t offsets[],
730 const LO numOffsets)
const;
739 virtual bool checkSizes (const ::Tpetra::SrcDistObject& source);
742 copyAndPermute (const ::Tpetra::SrcDistObject& source,
744 const Teuchos::ArrayView<const LO>& permuteToLIDs,
745 const Teuchos::ArrayView<const LO>& permuteFromLIDs);
748 packAndPrepare (const ::Tpetra::SrcDistObject& source,
749 const Teuchos::ArrayView<const LO>& exportLIDs,
750 Teuchos::Array<packet_type>& exports,
751 const Teuchos::ArrayView<size_t>& numPacketsPerLID,
752 size_t& constantNumPackets,
756 unpackAndCombine (
const Teuchos::ArrayView<const LO> &importLIDs,
757 const Teuchos::ArrayView<const packet_type> &imports,
758 const Teuchos::ArrayView<size_t> &numPacketsPerLID,
759 size_t constantNumPackets,
767 Teuchos::RCP<crs_graph_type> graphRCP_;
807 typename crs_graph_type::local_graph_type::row_map_type::HostMirror ptrHost_;
814 typename crs_graph_type::local_graph_type::entries_type::HostMirror indHost_;
821 typename Kokkos::DualView<impl_scalar_type*, device_type> val_;
844 Teuchos::RCP<Teuchos::RCP<BMV> > X_colMap_;
848 Teuchos::RCP<Teuchos::RCP<BMV> > Y_rowMap_;
857 Teuchos::RCP<Teuchos::RCP<typename crs_graph_type::import_type> > pointImporter_;
873 Teuchos::RCP<bool> localError_;
882 Teuchos::RCP<Teuchos::RCP<std::ostringstream> > errs_;
885 std::ostream& markLocalErrorAndGetStream ();
895 template<
class MemorySpace>
904 val_.template modify<typename MemorySpace::memory_space> ();
908 template<
class MemorySpace>
917 return val_.template need_sync<typename MemorySpace::memory_space> ();
921 template<
class MemorySpace>
930 val_.template sync<typename MemorySpace::memory_space> ();
947 template<
class MemorySpace>
948 auto getValues () -> decltype (val_.template view<typename MemorySpace::memory_space> ())
950 return val_.template view<typename MemorySpace::memory_space> ();
969 const Teuchos::ETransp mode,
1039 findRelOffsetOfColumnIndex (
const LO localRowIndex,
1040 const LO colIndexToFind,
1041 const LO hint = 0)
const;
1045 LO offsetPerBlock ()
const;
1048 getConstLocalBlockFromInput (
const impl_scalar_type* val,
const size_t pointOffset)
const;
1051 getNonConstLocalBlockFromInput (
impl_scalar_type* val,
const size_t pointOffset)
const;
1054 getConstLocalBlockFromAbsOffset (
const size_t absBlockOffset)
const;
1057 getNonConstLocalBlockFromAbsOffset (
const size_t absBlockOffset)
const;
1063 getConstLocalBlockFromRelOffset (
const LO lclMeshRow,
1064 const size_t relMeshOffset)
const;
1068 virtual Teuchos::RCP<const Teuchos::Comm<int> >
getComm()
const;
1071 virtual Teuchos::RCP<Node>
getNode()
const;
1199 const Teuchos::ArrayView<GO> &Indices,
1200 const Teuchos::ArrayView<Scalar> &Values,
1201 size_t& NumEntries)
const;
1229 Teuchos::ArrayView<const GO>& indices,
1230 Teuchos::ArrayView<const Scalar>& values)
const;
1254 virtual void leftScale (const ::Tpetra::Vector<Scalar, LO, GO, Node>& x);
1261 virtual void rightScale (const ::Tpetra::Vector<Scalar, LO, GO, Node>& x);
1271 virtual typename ::Tpetra::RowMatrix<Scalar, LO, GO, Node>::mag_type
1280 #endif // TPETRA_EXPERIMENTAL_BLOCKCRSMATRIX_DECL_HPP
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Print a description of this object to the given output stream.
virtual size_t getNodeNumCols() const
The number of columns needed to apply the forward operator on this node.
ESweepDirection
Sweep direction for Gauss-Seidel or Successive Over-Relaxation (SOR).
LO replaceLocalValues(const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
Replace values at the given (mesh, i.e., block) column indices, in the given (mesh,...
LO sumIntoLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like sumIntoLocalValues, but avoids computing row offsets.
virtual bool isFillComplete() const
Whether fillComplete() has been called.
virtual global_size_t getGlobalNumEntries() const
The global number of stored (structurally nonzero) entries.
LO absMaxLocalValues(const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
Like sumIntoLocalValues, but for the ABSMAX combine mode.
Kokkos::View< const impl_scalar_type **, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > const_little_block_type
The type used to access const matrix blocks.
virtual size_t getGlobalMaxNumRowEntries() const
The maximum number of entries across all rows/columns on all nodes.
std::string description() const
One-line description of this object.
Kokkos::View< const impl_scalar_type *, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > const_little_vec_type
The type used to access const vector blocks.
virtual Teuchos::RCP< const ::Tpetra::RowGraph< LO, GO, Node > > getGraph() const
Get the (mesh) graph.
void sync()
Sync the matrix's values to the given memory space.
bool need_sync() const
Whether the matrix's values need sync'ing to the given memory space.
Teuchos::RCP< const map_type > getRowMap() const
get the (mesh) map for the rows of this block matrix.
LO getLocalRowView(const LO localRowInd, const LO *&colInds, Scalar *&vals, LO &numInds) const
Get a view of the (mesh, i.e., block) row, using local (mesh, i.e., block) indices.
LO replaceLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like replaceLocalValues, but avoids computing row offsets.
::Kokkos::Compat::KokkosDeviceWrapperNode< execution_space > node_type
Default value of Node template parameter.
void apply(const mv_type &X, mv_type &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
For this matrix A, compute Y := beta * Y + alpha * Op(A) * X.
LO getBlockSize() const
The number of degrees of freedom per mesh point.
void modify()
Mark the matrix's values as modified in the given memory space.
auto getValues() -> decltype(val_.template view< typename MemorySpace::memory_space >())
Get the host or device View of the matrix's values (val_).
std::string errorMessages() const
The current stream of error messages.
virtual void getGlobalRowView(GO GlobalRow, Teuchos::ArrayView< const GO > &indices, Teuchos::ArrayView< const Scalar > &values) const
Get a constant, nonpersisting, globally indexed view of the given row of the matrix.
virtual bool TPETRA_DEPRECATED isUpperTriangular() const
Whether the matrix's graph is locally upper triangular.
Kokkos::View< impl_scalar_type *, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_vec_type
The type used to access nonconst vector blocks.
virtual void leftScale(const ::Tpetra::Vector< Scalar, LO, GO, Node > &x)
Scale the RowMatrix on the left with the given Vector x.
void getLocalRowCopy(LO LocalRow, const Teuchos::ArrayView< LO > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const
Not implemented.
virtual GO getIndexBase() const
The index base for global indices in this matrix.
size_t getNodeMaxNumRowEntries() const
Maximum number of entries in any row of the matrix, on this process.
char packet_type
Implementation detail; tells.
Sparse matrix whose entries are small dense square blocks, all of the same dimensions.
void localGaussSeidel(const BlockMultiVector< Scalar, LO, GO, Node > &Residual, BlockMultiVector< Scalar, LO, GO, Node > &Solution, const Kokkos::View< impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &D_inv, const Scalar &omega, const ESweepDirection direction) const
Local Gauss-Seidel solve, given a factorized diagonal.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to alpha.
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which this matrix is distributed.
Teuchos::RCP< const map_type > getColMap() const
get the (mesh) map for the columns of this block matrix.
mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the matrix.
Scalar scalar_type
The type of entries in the matrix (that is, of each entry in each block).
BlockCrsMatrix()
Default constructor: Makes an empty block matrix.
Kokkos::View< impl_scalar_type **, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_block_type
The type used to access nonconst matrix blocks.
Declaration of the Tpetra::CrsMatrix class.
typename BMV::impl_scalar_type impl_scalar_type
The implementation type of entries in the matrix.
virtual size_t TPETRA_DEPRECATED getNodeNumDiags() const
Number of diagonal entries in the matrix's graph, on the calling process.
Base class for distributed Tpetra objects that support data redistribution.
virtual void rightScale(const ::Tpetra::Vector< Scalar, LO, GO, Node > &x)
Scale the RowMatrix on the right with the given Vector x.
size_t getNodeNumRows() const
get the local number of block rows
Node node_type
The Node type.
virtual bool hasColMap() const
Whether this matrix has a well-defined column map.
Sets up and executes a communication plan for a Tpetra DistObject.
virtual typename ::Tpetra::RowMatrix< Scalar, LO, GO, Node >::mag_type getFrobeniusNorm() const
The Frobenius norm of the matrix.
::Tpetra::MultiVector< Scalar, LO, GO, node_type > mv_type
The implementation of MultiVector that this class uses.
virtual ~BlockCrsMatrix()
Destructor (declared virtual for memory safety).
One or more distributed dense vectors.
int local_ordinal_type
Default value of Scalar template parameter.
device_type::execution_space execution_space
The Kokkos execution space that this class uses.
virtual bool supportsRowViews() const
Whether this object implements getLocalRowView() and getGlobalRowView().
virtual global_size_t TPETRA_DEPRECATED getGlobalNumDiags() const
Number of diagonal entries in the matrix's graph, over all processes in the matrix's communicator.
void gaussSeidelCopy(MultiVector< Scalar, LO, GO, Node > &X, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &B, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &D, const Scalar &dampingFactor, const ESweepDirection direction, const int numSweeps, const bool zeroInitialGuess) const
Version of gaussSeidel(), with fewer requirements on X.
bool localError() const
Whether this object had an error on the calling process.
virtual size_t getNodeNumEntries() const
The local number of stored (structurally nonzero) entries.
virtual void getGlobalRowCopy(GO GlobalRow, const Teuchos::ArrayView< GO > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const
Get a copy of the given global row's entries.
virtual size_t getNumEntriesInGlobalRow(GO globalRow) const
The current number of entries on the calling process in the specified global row.
Teuchos::RCP< const map_type > getDomainMap() const
Get the (point) domain Map of this matrix.
GO global_ordinal_type
The type of global indices.
::Tpetra::Map< LO, GO, node_type > map_type
The implementation of Map that this class uses.
virtual bool isLocallyIndexed() const
Whether matrix indices are locally indexed.
virtual bool TPETRA_DEPRECATED isLowerTriangular() const
Whether the matrix's graph is locally lower triangular.
Node::device_type device_type
The Kokkos::Device specialization that this class uses.
MultiVector for multiple degrees of freedom per mesh point.
A parallel distribution of indices over processes.
LO absMaxLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like sumIntoLocalValuesByOffsets, but for the ABSMAX combine mode.
void getLocalDiagOffsets(const Kokkos::View< size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
Get offsets of the diagonal entries in the matrix.
void applyBlock(const BlockMultiVector< Scalar, LO, GO, Node > &X, BlockMultiVector< Scalar, LO, GO, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, const Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), const Scalar beta=Teuchos::ScalarTraits< Scalar >::zero())
Version of apply() that takes BlockMultiVector input and output.
void getLocalDiagCopy(const Kokkos::View< impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &diag, const Kokkos::View< const size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
Variant of getLocalDiagCopy() that uses precomputed offsets and puts diagonal blocks in a 3-D Kokkos:...
LO local_ordinal_type
The type of local indices.
void reorderedGaussSeidelCopy(::Tpetra::MultiVector< Scalar, LO, GO, Node > &X, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &B, const ::Tpetra::MultiVector< Scalar, LO, GO, Node > &D, const Teuchos::ArrayView< LO > &rowIndices, const Scalar &dampingFactor, const ESweepDirection direction, const int numSweeps, const bool zeroInitialGuess) const
Version of reorderedGaussSeidel(), with fewer requirements on X.
size_t global_size_t
Global size_t object.
bool hasTransposeApply() const
Whether it is valid to apply the transpose or conjugate transpose of this matrix.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
LO sumIntoLocalValues(const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
Sum into values at the given (mesh, i.e., block) column indices, in the given (mesh,...
global_size_t getGlobalNumRows() const
get the global number of block rows
Teuchos::RCP< const map_type > getRangeMap() const
Get the (point) range Map of this matrix.
virtual global_size_t getGlobalNumCols() const
The global number of columns of this matrix.
LO getLocalRowOffsets(const LO localRowInd, ptrdiff_t offsets[], const LO colInds[], const LO numColInds) const
Get relative offsets corresponding to the given rows, given by local row index.
device_type::memory_space memory_space
The Kokkos memory space that this class uses.
::Tpetra::CrsGraph< LO, GO, node_type > crs_graph_type
The implementation of CrsGraph that this class uses.
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
virtual bool isGloballyIndexed() const
Whether matrix indices are globally indexed.
size_t getNumEntriesInLocalRow(const LO localRowInd) const
Return the number of entries in the given row on the calling process.
virtual Teuchos::RCP< Node > getNode() const
The Kokkos Node instance.
CombineMode
Rule for combining data in an Import or Export.
A read-only, row-oriented interface to a sparse matrix.
A distributed dense vector.