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 //
4 // Tpetra: Templated Linear Algebra Services Package
5 // Copyright (2008) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 // @HEADER
41 
42 #ifndef TPETRA_ROWMATRIX_DECL_HPP
43 #define TPETRA_ROWMATRIX_DECL_HPP
44 
45 #include "Tpetra_ConfigDefs.hpp"
46 #include "Tpetra_RowMatrix_fwd.hpp"
47 #include "Tpetra_Vector_fwd.hpp"
48 #include "Tpetra_Operator.hpp"
49 #include "Tpetra_RowGraph_fwd.hpp"
50 #include "Tpetra_Packable.hpp"
51 #include "Tpetra_SrcDistObject.hpp"
52 #include "Teuchos_Describable.hpp"
53 #include "Kokkos_ArithTraits.hpp"
54 
55 namespace Tpetra {
81  template <class Scalar,
82  class LocalOrdinal,
83  class GlobalOrdinal,
84  class Node>
85  class RowMatrix :
86  virtual public Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
87  virtual public SrcDistObject,
88  public Packable<char, LocalOrdinal> {
89  public:
91 
92 
94  typedef Scalar scalar_type;
96  typedef LocalOrdinal local_ordinal_type;
98  typedef GlobalOrdinal global_ordinal_type;
100  typedef Node node_type;
101 
111  using impl_scalar_type = typename Kokkos::ArithTraits<Scalar>::val_type;
117  using mag_type = typename Kokkos::ArithTraits<Scalar>::mag_type;
118 
119  typedef typename
120  Kokkos::View<impl_scalar_type*, typename Node::device_type>::const_type
121  values_device_view_type;
122  typedef typename values_device_view_type::HostMirror::const_type
123  values_host_view_type;
124  typedef typename values_device_view_type::HostMirror
125  nonconst_values_host_view_type;
126 
127  typedef typename
128  Kokkos::View<LocalOrdinal *, typename Node::device_type>::const_type
129  local_inds_device_view_type;
130  typedef typename local_inds_device_view_type::HostMirror::const_type
131  local_inds_host_view_type;
132  typedef typename local_inds_device_view_type::HostMirror
133  nonconst_local_inds_host_view_type;
134 
135  typedef typename
136  Kokkos::View<GlobalOrdinal *, typename Node::device_type>::const_type
137  global_inds_device_view_type;
138  typedef typename global_inds_device_view_type::HostMirror::const_type
139  global_inds_host_view_type;
140  typedef typename global_inds_device_view_type::HostMirror
141  nonconst_global_inds_host_view_type;
142 
143 
144  typedef typename
145  Kokkos::View<const size_t*, typename Node::device_type>::const_type
146  row_ptrs_device_view_type;
147  typedef typename row_ptrs_device_view_type::HostMirror::const_type
148  row_ptrs_host_view_type;
149 
150 
151 
153 
155 
157  virtual ~RowMatrix();
158 
160 
162 
164  virtual Teuchos::RCP<const Teuchos::Comm<int> > getComm() const = 0;
165 
166 
168  virtual Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > getRowMap() const = 0;
169 
171  virtual Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > getColMap() const = 0;
172 
174  virtual Teuchos::RCP<const RowGraph<LocalOrdinal,GlobalOrdinal,Node> > getGraph() const = 0;
175 
177  virtual global_size_t getGlobalNumRows() const = 0;
178 
180  virtual global_size_t getGlobalNumCols() const = 0;
181 
183  virtual size_t getLocalNumRows() const = 0;
184 
190  virtual size_t getLocalNumCols() const = 0;
191 
193  virtual GlobalOrdinal getIndexBase() const = 0;
194 
196  virtual global_size_t getGlobalNumEntries() const = 0;
197 
199  virtual size_t getLocalNumEntries() const = 0;
200 
210  virtual size_t getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const = 0;
211 
221  virtual size_t getNumEntriesInLocalRow (LocalOrdinal localRow) const = 0;
222 
231  virtual size_t getGlobalMaxNumRowEntries () const = 0;
232 
234  virtual local_ordinal_type getBlockSize () const = 0;
235 
244  virtual size_t getLocalMaxNumRowEntries () const = 0;
245 
247  virtual bool hasColMap () const = 0;
248 
258  virtual bool isLocallyIndexed() const = 0;
259 
269  virtual bool isGloballyIndexed() const = 0;
270 
272  virtual bool isFillComplete() const = 0;
273 
275  virtual bool supportsRowViews() const = 0;
276 
277 
279 
281 
302  virtual void
303  getGlobalRowCopy (GlobalOrdinal GlobalRow,
304  nonconst_global_inds_host_view_type &Indices,
305  nonconst_values_host_view_type &Values,
306  size_t& NumEntries) const = 0;
307 
328  virtual void
329  getLocalRowCopy (LocalOrdinal LocalRow,
330  nonconst_local_inds_host_view_type &Indices,
331  nonconst_values_host_view_type &Values,
332  size_t& NumEntries) const = 0;
333 
358  virtual void
359  getGlobalRowView (GlobalOrdinal GlobalRow,
360  global_inds_host_view_type &indices,
361  values_host_view_type &values) const = 0;
362 
387  virtual void
388  getLocalRowView (LocalOrdinal LocalRow,
389  local_inds_host_view_type & indices,
390  values_host_view_type & values) const = 0;
391 
404 
406 
408 
414 
420 
429  virtual mag_type getFrobeniusNorm() const = 0;
430 
482  virtual Teuchos::RCP<RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
483  add (const Scalar& alpha,
485  const Scalar& beta,
486  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& domainMap = Teuchos::null,
487  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& rangeMap = Teuchos::null,
488  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const;
490 
492  private:
493  bool
494  packRow (char* const numEntOut,
495  char* const valOut,
496  char* const indOut,
497  const size_t numEnt,
498  const LocalOrdinal lclRow) const;
499 
500  // TODO (mfh 25 Jan 2015) Could just make this "protected" and let
501  // CrsMatrix use it, since it's exactly the same there.
502  void
503  allocatePackSpace (Teuchos::Array<char>& exports,
504  size_t& totalNumEntries,
505  const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs) const;
506 
511  void
512  packImpl (const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs,
513  Teuchos::Array<char>& exports,
514  const Teuchos::ArrayView<size_t>& numPacketsPerLID,
515  size_t& constantNumPackets) const;
516 
517 
518  public:
527  virtual void
528  pack (const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs,
529  Teuchos::Array<char>& exports,
530  const Teuchos::ArrayView<size_t>& numPacketsPerLID,
531  size_t& constantNumPackets) const;
533  }; // class RowMatrix
534 } // namespace Tpetra
535 
536 #endif // TPETRA_ROWMATRIX_DECL_HPP
537 
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.