Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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 
107  using mag_type = typename Kokkos::ArithTraits<Scalar>::mag_type;
108 
110 
112 
114  virtual ~RowMatrix();
115 
117 
119 
121  virtual Teuchos::RCP<const Teuchos::Comm<int> > getComm() const = 0;
122 
123 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
124  virtual TPETRA_DEPRECATED Teuchos::RCP<Node> getNode() const = 0;
126 #endif // TPETRA_ENABLE_DEPRECATED_CODE
127 
129  virtual Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > getRowMap() const = 0;
130 
132  virtual Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > getColMap() const = 0;
133 
135  virtual Teuchos::RCP<const RowGraph<LocalOrdinal,GlobalOrdinal,Node> > getGraph() const = 0;
136 
138  virtual global_size_t getGlobalNumRows() const = 0;
139 
141  virtual global_size_t getGlobalNumCols() const = 0;
142 
144  virtual size_t getNodeNumRows() const = 0;
145 
151  virtual size_t getNodeNumCols() const = 0;
152 
154  virtual GlobalOrdinal getIndexBase() const = 0;
155 
157  virtual global_size_t getGlobalNumEntries() const = 0;
158 
160  virtual size_t getNodeNumEntries() const = 0;
161 
171  virtual size_t getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const = 0;
172 
182  virtual size_t getNumEntriesInLocalRow (LocalOrdinal localRow) const = 0;
183 
192  virtual size_t getGlobalMaxNumRowEntries () const = 0;
193 
202  virtual size_t getNodeMaxNumRowEntries () const = 0;
203 
205  virtual bool hasColMap () const = 0;
206 
216  virtual bool isLocallyIndexed() const = 0;
217 
227  virtual bool isGloballyIndexed() const = 0;
228 
230  virtual bool isFillComplete() const = 0;
231 
233  virtual bool supportsRowViews() const = 0;
234 
235 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
236  virtual global_size_t TPETRA_DEPRECATED getGlobalNumDiags () const = 0;
245 
254  virtual size_t TPETRA_DEPRECATED getNodeNumDiags () const = 0;
255 
266  virtual bool TPETRA_DEPRECATED isLowerTriangular () const = 0;
267 
278  virtual bool TPETRA_DEPRECATED isUpperTriangular () const = 0;
279 #endif // TPETRA_ENABLE_DEPRECATED_CODE
280 
282 
284 
305  virtual void
306  getGlobalRowCopy (GlobalOrdinal GlobalRow,
307  const Teuchos::ArrayView<GlobalOrdinal> &Indices,
308  const Teuchos::ArrayView<Scalar> &Values,
309  size_t &NumEntries) const = 0;
310 
331  virtual void
332  getLocalRowCopy (LocalOrdinal LocalRow,
333  const Teuchos::ArrayView<LocalOrdinal> &Indices,
334  const Teuchos::ArrayView<Scalar> &Values,
335  size_t &NumEntries) const = 0;
336 
361  virtual void
362  getGlobalRowView (GlobalOrdinal GlobalRow,
363  Teuchos::ArrayView<const GlobalOrdinal> &indices,
364  Teuchos::ArrayView<const Scalar> &values) const = 0;
365 
390  virtual void
391  getLocalRowView (LocalOrdinal LocalRow,
392  Teuchos::ArrayView<const LocalOrdinal>& indices,
393  Teuchos::ArrayView<const Scalar>& values) const = 0;
394 
421  virtual LocalOrdinal
422  getLocalRowViewRaw (const LocalOrdinal lclRow,
423  LocalOrdinal& numEnt,
424  const LocalOrdinal*& lclColInds,
425  const Scalar*& vals) const;
426 
439 
441 
443 
449 
455 
464  virtual mag_type getFrobeniusNorm() const = 0;
465 
517  virtual Teuchos::RCP<RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
518  add (const Scalar& alpha,
520  const Scalar& beta,
521  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& domainMap = Teuchos::null,
522  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >& rangeMap = Teuchos::null,
523  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) const;
525 
527  private:
528  bool
529  packRow (char* const numEntOut,
530  char* const valOut,
531  char* const indOut,
532  const size_t numEnt,
533  const LocalOrdinal lclRow) const;
534 
535  // TODO (mfh 25 Jan 2015) Could just make this "protected" and let
536  // CrsMatrix use it, since it's exactly the same there.
537  void
538  allocatePackSpace (Teuchos::Array<char>& exports,
539  size_t& totalNumEntries,
540  const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs) const;
541 
546  void
547  packImpl (const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs,
548  Teuchos::Array<char>& exports,
549  const Teuchos::ArrayView<size_t>& numPacketsPerLID,
550  size_t& constantNumPackets,
551  Distributor& distor) const;
552 
553 
554  public:
563  virtual void
564  pack (const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs,
565  Teuchos::Array<char>& exports,
566  const Teuchos::ArrayView<size_t>& numPacketsPerLID,
567  size_t& constantNumPackets,
568  Distributor& distor) const;
570  }; // class RowMatrix
571 } // namespace Tpetra
572 
573 #endif // TPETRA_ROWMATRIX_DECL_HPP
574 
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 size_t getNodeMaxNumRowEntries() const =0
Maximum number of entries in any row of the matrix, on this process.
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 size_t getNodeNumRows() const =0
The number of rows owned by the calling process.
virtual void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Scale the matrix on the right with the given Vector.
virtual size_t getNodeNumCols() const =0
The number of columns needed to apply the forward operator on this node.
virtual void getLocalRowCopy(LocalOrdinal LocalRow, const Teuchos::ArrayView< LocalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const =0
Get a copy of the given local row&#39;s entries.
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 getLocalRowView(LocalOrdinal LocalRow, Teuchos::ArrayView< const LocalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const =0
Get a constant, nonpersisting, locally indexed view of the given row of the matrix.
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 LocalOrdinal getLocalRowViewRaw(const LocalOrdinal lclRow, LocalOrdinal &numEnt, const LocalOrdinal *&lclColInds, const Scalar *&vals) const
Get a constant, nonpersisting, locally indexed view of the given row of the matrix, using &quot;raw&quot; pointers instead of Teuchos::ArrayView.
virtual ~RowMatrix()
Destructor (virtual for memory safety of derived classes).
Forward declaration of Tpetra::RowGraph.
size_t global_size_t
Global size_t object.
virtual void getGlobalRowView(GlobalOrdinal GlobalRow, Teuchos::ArrayView< const GlobalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const =0
Get a constant, nonpersisting, globally indexed view of the given row of the matrix.
virtual void getGlobalRowCopy(GlobalOrdinal GlobalRow, const Teuchos::ArrayView< GlobalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const =0
Get a copy of the given global row&#39;s entries.
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).
typename Kokkos::ArithTraits< Scalar >::mag_type mag_type
Type of a norm result.
GlobalOrdinal global_ordinal_type
The type of global indices.
virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const =0
The current number of entries on the calling process in the specified global row. ...
Sets up and executes a communication plan for a Tpetra DistObject.
virtual void pack(const Teuchos::ArrayView< const LocalOrdinal > &exportLIDs, Teuchos::Array< char > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets, Distributor &distor) const
Pack this object&#39;s data for an Import or Export.
CrsGraphType::global_ordinal_type getGlobalNumDiags(const CrsGraphType &G)
Number of populated diagonal entries in the given sparse graph, over all processes in the graph&#39;s (MP...
Abstract base class for objects that can be the source of an Import or Export operation.
Forward declaration of Tpetra::Vector.
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 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 size_t getNodeNumEntries() const =0
The local number of stored (structurally nonzero) entries.
virtual bool isGloballyIndexed() const =0
Whether matrix indices are globally indexed.
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 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 mag_type getFrobeniusNorm() const =0
The Frobenius norm of the matrix.