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 #include "Kokkos_ArithTraits.hpp"
22 
23 namespace Tpetra {
49  template <class Scalar,
50  class LocalOrdinal,
51  class GlobalOrdinal,
52  class Node>
53  class RowMatrix :
54  virtual public Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
55  virtual public SrcDistObject,
56  public Packable<char, LocalOrdinal> {
57  public:
59 
60 
62  typedef Scalar scalar_type;
64  typedef LocalOrdinal local_ordinal_type;
66  typedef GlobalOrdinal global_ordinal_type;
68  typedef Node node_type;
69 
79  using impl_scalar_type = typename Kokkos::ArithTraits<Scalar>::val_type;
85  using mag_type = typename Kokkos::ArithTraits<Scalar>::mag_type;
86 
87  typedef typename
88  Kokkos::View<impl_scalar_type*, typename Node::device_type>::const_type
89  values_device_view_type;
90  typedef typename values_device_view_type::HostMirror::const_type
91  values_host_view_type;
92  typedef typename values_device_view_type::HostMirror
93  nonconst_values_host_view_type;
94 
95  typedef typename
96  Kokkos::View<LocalOrdinal *, typename Node::device_type>::const_type
97  local_inds_device_view_type;
98  typedef typename local_inds_device_view_type::HostMirror::const_type
99  local_inds_host_view_type;
100  typedef typename local_inds_device_view_type::HostMirror
101  nonconst_local_inds_host_view_type;
102 
103  typedef typename
104  Kokkos::View<GlobalOrdinal *, typename Node::device_type>::const_type
105  global_inds_device_view_type;
106  typedef typename global_inds_device_view_type::HostMirror::const_type
107  global_inds_host_view_type;
108  typedef typename global_inds_device_view_type::HostMirror
109  nonconst_global_inds_host_view_type;
110 
111 
112  typedef typename
113  Kokkos::View<const size_t*, typename Node::device_type>::const_type
114  row_ptrs_device_view_type;
115  typedef typename row_ptrs_device_view_type::HostMirror::const_type
116  row_ptrs_host_view_type;
117 
118 
119 
121 
123 
125  virtual ~RowMatrix();
126 
128 
130 
132  virtual Teuchos::RCP<const Teuchos::Comm<int> > getComm() const = 0;
133 
134 
136  virtual Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > getRowMap() const = 0;
137 
139  virtual Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > getColMap() const = 0;
140 
142  virtual Teuchos::RCP<const RowGraph<LocalOrdinal,GlobalOrdinal,Node> > getGraph() const = 0;
143 
145  virtual global_size_t getGlobalNumRows() const = 0;
146 
148  virtual global_size_t getGlobalNumCols() const = 0;
149 
151  virtual size_t getLocalNumRows() const = 0;
152 
158  virtual size_t getLocalNumCols() const = 0;
159 
161  virtual GlobalOrdinal getIndexBase() const = 0;
162 
164  virtual global_size_t getGlobalNumEntries() const = 0;
165 
167  virtual size_t getLocalNumEntries() const = 0;
168 
178  virtual size_t getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const = 0;
179 
189  virtual size_t getNumEntriesInLocalRow (LocalOrdinal localRow) const = 0;
190 
199  virtual size_t getGlobalMaxNumRowEntries () const = 0;
200 
202  virtual local_ordinal_type getBlockSize () const = 0;
203 
212  virtual size_t getLocalMaxNumRowEntries () const = 0;
213 
215  virtual bool hasColMap () const = 0;
216 
226  virtual bool isLocallyIndexed() const = 0;
227 
237  virtual bool isGloballyIndexed() const = 0;
238 
240  virtual bool isFillComplete() const = 0;
241 
243  virtual bool supportsRowViews() const = 0;
244 
245 
247 
249 
270  virtual void
271  getGlobalRowCopy (GlobalOrdinal GlobalRow,
272  nonconst_global_inds_host_view_type &Indices,
273  nonconst_values_host_view_type &Values,
274  size_t& NumEntries) const = 0;
275 
296  virtual void
297  getLocalRowCopy (LocalOrdinal LocalRow,
298  nonconst_local_inds_host_view_type &Indices,
299  nonconst_values_host_view_type &Values,
300  size_t& NumEntries) const = 0;
301 
326  virtual void
327  getGlobalRowView (GlobalOrdinal GlobalRow,
328  global_inds_host_view_type &indices,
329  values_host_view_type &values) const = 0;
330 
355  virtual void
356  getLocalRowView (LocalOrdinal LocalRow,
357  local_inds_host_view_type & indices,
358  values_host_view_type & values) const = 0;
359 
372 
374 
376 
382 
388 
397  virtual mag_type getFrobeniusNorm() const = 0;
398 
450  virtual Teuchos::RCP<RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
451  add (const Scalar& alpha,
453  const Scalar& beta,
454  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& domainMap = Teuchos::null,
455  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& rangeMap = Teuchos::null,
456  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const;
458 
460  private:
461  bool
462  packRow (char* const numEntOut,
463  char* const valOut,
464  char* const indOut,
465  const size_t numEnt,
466  const LocalOrdinal lclRow) const;
467 
468  // TODO (mfh 25 Jan 2015) Could just make this "protected" and let
469  // CrsMatrix use it, since it's exactly the same there.
470  void
471  allocatePackSpace (Teuchos::Array<char>& exports,
472  size_t& totalNumEntries,
473  const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs) const;
474 
479  void
480  packImpl (const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs,
481  Teuchos::Array<char>& exports,
482  const Teuchos::ArrayView<size_t>& numPacketsPerLID,
483  size_t& constantNumPackets) const;
484 
485 
486  public:
495  virtual void
496  pack (const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs,
497  Teuchos::Array<char>& exports,
498  const Teuchos::ArrayView<size_t>& numPacketsPerLID,
499  size_t& constantNumPackets) const;
501  }; // class RowMatrix
502 } // namespace Tpetra
503 
504 #endif // TPETRA_ROWMATRIX_DECL_HPP
505 
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const =0
The communicator over which this matrix is distributed.
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const =0
The Map that describes the distribution of columns over processes.
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.
virtual Teuchos::RCP< const RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const =0
The RowGraph associated with this matrix.
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.
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 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 Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const =0
The Map that describes the distribution of rows over processes.
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.