Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_RowMatrix_decl.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Tpetra: Templated Linear Algebra Services Package
4 //
5 // Copyright 2008 NTESS and the Tpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef TPETRA_ROWMATRIX_DECL_HPP
11 #define TPETRA_ROWMATRIX_DECL_HPP
12 
13 #include "Tpetra_ConfigDefs.hpp"
14 #include "Tpetra_RowMatrix_fwd.hpp"
15 #include "Tpetra_Vector_fwd.hpp"
16 #include "Tpetra_Operator.hpp"
17 #include "Tpetra_RowGraph_fwd.hpp"
18 #include "Tpetra_Packable.hpp"
19 #include "Tpetra_SrcDistObject.hpp"
20 #include "Teuchos_Describable.hpp"
21 #if KOKKOS_VERSION >= 40799
22 #include "KokkosKernels_ArithTraits.hpp"
23 #else
24 #include "Kokkos_ArithTraits.hpp"
25 #endif
26 
27 namespace Tpetra {
53 template <class Scalar,
54  class LocalOrdinal,
55  class GlobalOrdinal,
56  class Node>
57 class RowMatrix : virtual public Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
58  virtual public SrcDistObject,
59  public Packable<char, LocalOrdinal> {
60  public:
62 
63 
65  typedef Scalar scalar_type;
67  typedef LocalOrdinal local_ordinal_type;
69  typedef GlobalOrdinal global_ordinal_type;
71  typedef Node node_type;
72 
82 #if KOKKOS_VERSION >= 40799
83  using impl_scalar_type = typename KokkosKernels::ArithTraits<Scalar>::val_type;
84 #else
85  using impl_scalar_type = typename Kokkos::ArithTraits<Scalar>::val_type;
86 #endif
87 #if KOKKOS_VERSION >= 40799
93  using mag_type = typename KokkosKernels::ArithTraits<Scalar>::mag_type;
94 #else
95  using mag_type = typename Kokkos::ArithTraits<Scalar>::mag_type;
96 #endif
97 
98  typedef typename Kokkos::View<impl_scalar_type*, typename Node::device_type>::const_type
99  values_device_view_type;
100  typedef typename values_device_view_type::host_mirror_type::const_type
101  values_host_view_type;
102  typedef typename values_device_view_type::host_mirror_type
103  nonconst_values_host_view_type;
104 
105  typedef typename Kokkos::View<LocalOrdinal*, typename Node::device_type>::const_type
106  local_inds_device_view_type;
107  typedef typename local_inds_device_view_type::host_mirror_type::const_type
108  local_inds_host_view_type;
109  typedef typename local_inds_device_view_type::host_mirror_type
110  nonconst_local_inds_host_view_type;
111 
112  typedef typename Kokkos::View<GlobalOrdinal*, typename Node::device_type>::const_type
113  global_inds_device_view_type;
114  typedef typename global_inds_device_view_type::host_mirror_type::const_type
115  global_inds_host_view_type;
116  typedef typename global_inds_device_view_type::host_mirror_type
117  nonconst_global_inds_host_view_type;
118 
119  typedef typename Kokkos::View<const size_t*, typename Node::device_type>::const_type
120  row_ptrs_device_view_type;
121  typedef typename row_ptrs_device_view_type::host_mirror_type::const_type
122  row_ptrs_host_view_type;
123 
125 
127 
129  virtual ~RowMatrix();
130 
132 
134 
136  virtual Teuchos::RCP<const Teuchos::Comm<int> > getComm() const = 0;
137 
139  virtual Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > getRowMap() const = 0;
140 
142  virtual Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > getColMap() const = 0;
143 
145  virtual Teuchos::RCP<const RowGraph<LocalOrdinal, GlobalOrdinal, Node> > getGraph() const = 0;
146 
148  virtual global_size_t getGlobalNumRows() const = 0;
149 
151  virtual global_size_t getGlobalNumCols() const = 0;
152 
154  virtual size_t getLocalNumRows() const = 0;
155 
161  virtual size_t getLocalNumCols() const = 0;
162 
164  virtual GlobalOrdinal getIndexBase() const = 0;
165 
167  virtual global_size_t getGlobalNumEntries() const = 0;
168 
170  virtual size_t getLocalNumEntries() const = 0;
171 
181  virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const = 0;
182 
192  virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const = 0;
193 
202  virtual size_t getGlobalMaxNumRowEntries() const = 0;
203 
205  virtual local_ordinal_type getBlockSize() const = 0;
206 
215  virtual size_t getLocalMaxNumRowEntries() const = 0;
216 
218  virtual bool hasColMap() const = 0;
219 
229  virtual bool isLocallyIndexed() const = 0;
230 
240  virtual bool isGloballyIndexed() const = 0;
241 
243  virtual bool isFillComplete() const = 0;
244 
246  virtual bool supportsRowViews() const = 0;
247 
249 
251 
272  virtual void
273  getGlobalRowCopy(GlobalOrdinal GlobalRow,
274  nonconst_global_inds_host_view_type& Indices,
275  nonconst_values_host_view_type& Values,
276  size_t& NumEntries) const = 0;
277 
298  virtual void
299  getLocalRowCopy(LocalOrdinal LocalRow,
300  nonconst_local_inds_host_view_type& Indices,
301  nonconst_values_host_view_type& Values,
302  size_t& NumEntries) const = 0;
303 
328  virtual void
329  getGlobalRowView(GlobalOrdinal GlobalRow,
330  global_inds_host_view_type& indices,
331  values_host_view_type& values) const = 0;
332 
357  virtual void
358  getLocalRowView(LocalOrdinal LocalRow,
359  local_inds_host_view_type& indices,
360  values_host_view_type& values) const = 0;
361 
374 
376 
378 
384 
390 
399  virtual mag_type getFrobeniusNorm() const = 0;
400 
452  virtual Teuchos::RCP<RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
453  add(const Scalar& alpha,
455  const Scalar& beta,
456  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& domainMap = Teuchos::null,
457  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& rangeMap = Teuchos::null,
458  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const;
460 
462  private:
463  bool
464  packRow(char* const numEntOut,
465  char* const valOut,
466  char* const indOut,
467  const size_t numEnt,
468  const LocalOrdinal lclRow) const;
469 
470  // TODO (mfh 25 Jan 2015) Could just make this "protected" and let
471  // CrsMatrix use it, since it's exactly the same there.
472  void
473  allocatePackSpace(Teuchos::Array<char>& exports,
474  size_t& totalNumEntries,
475  const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs) const;
476 
481  void
482  packImpl(const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs,
483  Teuchos::Array<char>& exports,
484  const Teuchos::ArrayView<size_t>& numPacketsPerLID,
485  size_t& constantNumPackets) const;
486 
487  public:
496  virtual void
497  pack(const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs,
498  Teuchos::Array<char>& exports,
499  const Teuchos::ArrayView<size_t>& numPacketsPerLID,
500  size_t& constantNumPackets) const;
502 }; // class RowMatrix
503 } // namespace Tpetra
504 
505 #endif // TPETRA_ROWMATRIX_DECL_HPP
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const =0
The communicator over which this matrix is distributed.
virtual bool isLocallyIndexed() const =0
Whether matrix indices are locally indexed.
virtual global_size_t getGlobalNumCols() const =0
The global number of columns of this matrix.
virtual GlobalOrdinal getIndexBase() const =0
The index base for global indices in this matrix.
Node node_type
The Kokkos Node type.
virtual void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Scale the matrix on the right with the given Vector.
Forward declaration of Tpetra::RowMatrix.
virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const =0
The current number of entries on the calling process in the specified local row.
virtual void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Scale the matrix on the left with the given Vector.
Abstract base class for sources of an Import or Export.
LocalOrdinal local_ordinal_type
The type of local indices.
virtual ~RowMatrix()
Destructor (virtual for memory safety of derived classes).
Forward declaration of Tpetra::RowGraph.
typename Kokkos::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
size_t global_size_t
Global size_t object.
virtual size_t getGlobalMaxNumRowEntries() const =0
Maximum number of entries in any row of the matrix, over all processes.
Abstract interface for operators (e.g., matrices and preconditioners).
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&#39;s entries.
typename Kokkos::ArithTraits< Scalar >::mag_type mag_type
Type of a norm result.
GlobalOrdinal global_ordinal_type
The type of global indices.
virtual local_ordinal_type getBlockSize() const =0
The number of degrees of freedom per mesh point.
virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const =0
The current number of entries on the calling process in the specified global row. ...
Abstract base class for objects that can be the source of an Import or Export operation.
virtual size_t getLocalNumCols() const =0
The number of columns needed to apply the forward operator on this node.
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const =0
The Map that describes the distribution of columns over processes.
Forward declaration of Tpetra::Vector.
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&#39;s data for an Import or Export.
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const =0
The Map that describes the distribution of rows over processes.
virtual void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const =0
Get a copy of the diagonal entries, distributed by the row Map.
virtual bool isFillComplete() const =0
Whether fillComplete() has been called.
virtual global_size_t getGlobalNumEntries() const =0
The global number of stored (structurally nonzero) entries.
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=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) const
Return a new RowMatrix which is the result of beta*this + alpha*A.
Scalar scalar_type
The type of the entries in the sparse matrix.
A parallel distribution of indices over processes.
virtual size_t getLocalMaxNumRowEntries() const =0
Maximum number of entries in any row of the matrix, on this process.
virtual bool supportsRowViews() const =0
Whether this object implements getLocalRowView() and getGlobalRowView().
A read-only, row-oriented interface to a sparse matrix.
A distributed dense vector.
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.
virtual bool isGloballyIndexed() const =0
Whether matrix indices are globally indexed.
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&#39;s entries.
Declaration of Tpetra::Packable.
Abstract base class for objects that can be the source of an Import or Export operation, and that also know how to pack their data to send to the target object.
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.
LO local_ordinal_type
The local index type.
virtual size_t getLocalNumEntries() const =0
The local number of stored (structurally nonzero) entries.
virtual Teuchos::RCP< const RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const =0
The RowGraph associated with this matrix.
virtual bool hasColMap() const =0
Whether this matrix has a well-defined column Map.
virtual global_size_t getGlobalNumRows() const =0
The global number of rows of this matrix.
virtual size_t getLocalNumRows() const =0
The number of rows owned by the calling process.
virtual mag_type getFrobeniusNorm() const =0
The Frobenius norm of the matrix.