Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_CrsMatrix.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Xpetra: A linear algebra interface package
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef XPETRA_CRSMATRIX_HPP
47 #define XPETRA_CRSMATRIX_HPP
48 
49 /* this file is automatically generated - do not edit (see script/interfaces.py) */
50 
51 #include <Tpetra_KokkosCompat_DefaultNode.hpp>
52 #include "Xpetra_ConfigDefs.hpp"
53 #include "Xpetra_RowMatrix.hpp"
54 #include "Xpetra_DistObject.hpp"
55 #include "Xpetra_CrsGraph.hpp"
56 #include "Xpetra_Vector.hpp"
57 
58 #ifdef HAVE_XPETRA_TPETRA
59 #include <Kokkos_StaticCrsGraph.hpp>
60 #include <KokkosSparse_CrsMatrix.hpp>
61 #endif
62 
63 namespace Xpetra {
64 
65 template <class Scalar,
66  class LocalOrdinal,
67  class GlobalOrdinal,
68  class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
69 class CrsMatrix
70  : public RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
71  public DistObject<char, LocalOrdinal, GlobalOrdinal, Node> {
72  public:
73  typedef Scalar scalar_type;
74  typedef LocalOrdinal local_ordinal_type;
75  typedef GlobalOrdinal global_ordinal_type;
76  typedef Node node_type;
77 
79 
80 
82  virtual ~CrsMatrix() {}
83 
85 
87 
89  virtual void
90  insertGlobalValues(GlobalOrdinal globalRow,
91  const ArrayView<const GlobalOrdinal> &cols,
92  const ArrayView<const Scalar> &vals) = 0;
93 
95  virtual void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals) = 0;
96 
98  virtual void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals) = 0;
99 
101  virtual void replaceLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals) = 0;
102 
104  virtual void setAllToScalar(const Scalar &alpha) = 0;
105 
107  virtual void scale(const Scalar &alpha) = 0;
108 
119  virtual void allocateAllValues(size_t numNonZeros, ArrayRCP<size_t> &rowptr, ArrayRCP<LocalOrdinal> &colind, ArrayRCP<Scalar> &values) = 0;
120 
122  virtual void setAllValues(const ArrayRCP<size_t> &rowptr, const ArrayRCP<LocalOrdinal> &colind, const ArrayRCP<Scalar> &values) = 0;
123 
125  virtual void getAllValues(ArrayRCP<const size_t> &rowptr, ArrayRCP<const LocalOrdinal> &colind, ArrayRCP<const Scalar> &values) const = 0;
126 
128  virtual void getAllValues(ArrayRCP<Scalar> &values) = 0;
129 
131 
133 
134 
136  virtual void resumeFill(const RCP<ParameterList> &params = null) = 0;
137 
139  virtual void fillComplete(const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &domainMap, const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &rangeMap, const RCP<ParameterList> &params = null) = 0;
140 
142  virtual void fillComplete(const RCP<ParameterList> &params = null) = 0;
143 
145  virtual void replaceDomainMapAndImporter(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &newDomainMap, Teuchos::RCP<const Import<LocalOrdinal, GlobalOrdinal, Node> > &newImporter) = 0;
146 
148  virtual void expertStaticFillComplete(const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &domainMap,
149  const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &rangeMap,
150  const RCP<const Import<LocalOrdinal, GlobalOrdinal, Node> > &importer = Teuchos::null,
151  const RCP<const Export<LocalOrdinal, GlobalOrdinal, Node> > &exporter = Teuchos::null,
152  const RCP<ParameterList> &params = Teuchos::null) = 0;
154 
156 
157 
159  virtual const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > getRowMap() const = 0;
160 
162  virtual const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > getColMap() const = 0;
163 
165  virtual RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node> > getCrsGraph() const = 0;
166 
168  virtual global_size_t getGlobalNumRows() const = 0;
169 
171  virtual global_size_t getGlobalNumCols() const = 0;
172 
174  virtual size_t getLocalNumRows() const = 0;
175 
177  virtual global_size_t getGlobalNumEntries() const = 0;
178 
180  virtual size_t getLocalNumEntries() const = 0;
181 
183  virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const = 0;
184 
186  virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const = 0;
187 
189  virtual size_t getGlobalMaxNumRowEntries() const = 0;
190 
192  virtual size_t getLocalMaxNumRowEntries() const = 0;
193 
195  virtual bool isLocallyIndexed() const = 0;
196 
198  virtual bool isGloballyIndexed() const = 0;
199 
201  virtual bool isFillComplete() const = 0;
202 
204  virtual bool isFillActive() const = 0;
205 
207  virtual typename ScalarTraits<Scalar>::magnitudeType getFrobeniusNorm() const = 0;
208 
210  virtual bool supportsRowViews() const = 0;
211 
213  virtual void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) const = 0;
214 
216  virtual void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView<GlobalOrdinal> &indices, const ArrayView<Scalar> &values, size_t &numEntries) const = 0;
217 
219  virtual void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) const = 0;
220 
223 
225  virtual void getLocalDiagOffsets(Teuchos::ArrayRCP<size_t> &offsets) const = 0;
226 
228  virtual void getLocalDiagCopy(Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &diag, const Teuchos::ArrayView<const size_t> &offsets) const = 0;
229 
231  virtual void getLocalDiagCopy(Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &diag, const Kokkos::View<const size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged> &offsets) const = 0;
232 
235 
238 
241 
242  virtual void removeEmptyProcessesInPlace(const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &newMap) = 0;
243 
245  virtual bool haveGlobalConstants() const = 0;
246 
248 
250 
251 
264  virtual void apply(const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &X, MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &Y, Teuchos::ETransp mode = Teuchos::NO_TRANS, Scalar alpha = ScalarTraits<Scalar>::one(), Scalar beta = ScalarTraits<Scalar>::zero()) const = 0;
265 
283  virtual void apply(const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &X, MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP<Import<LocalOrdinal, GlobalOrdinal, Node> > &regionInterfaceImporter, const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const = 0;
284 
286  virtual const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > getDomainMap() const = 0;
287 
289  virtual const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > getRangeMap() const = 0;
290 
292 
294 
295 
297  virtual std::string description() const = 0;
298 
300  virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default) const = 0;
301 
303 
305 
306  virtual void setObjectLabel(const std::string &objectLabel) = 0;
308 
310 
311 #ifdef HAVE_XPETRA_TPETRA
312  using impl_scalar_type = typename Kokkos::ArithTraits<Scalar>::val_type;
313  using execution_space = typename node_type::device_type;
314 
315  // that is the local_graph_type in Tpetra::CrsGraph...
316  using local_graph_type = Kokkos::StaticCrsGraph<LocalOrdinal,
317  Kokkos::LayoutLeft,
319  void,
320  size_t>;
324  using local_matrix_type = KokkosSparse::CrsMatrix<impl_scalar_type, LocalOrdinal, execution_space, void,
325  typename local_graph_type::size_type>;
326 
327  virtual local_matrix_type getLocalMatrixDevice() const = 0;
328  virtual typename local_matrix_type::HostMirror getLocalMatrixHost() const = 0;
329 
330  virtual void setAllValues(const typename local_matrix_type::row_map_type &ptr,
331  const typename local_graph_type::entries_type::non_const_type &ind,
332  const typename local_matrix_type::values_type &val) = 0;
333 #else
334 #ifdef __GNUC__
335 #warning "Xpetra Kokkos interface for CrsMatrix is enabled (HAVE_XPETRA_KOKKOS_REFACTOR) but Tpetra is disabled. The Kokkos interface needs Tpetra to be enabled, too."
336 #endif
337 #endif
338 
340 
341  // Adding these functions by hand, as they're in the skip list.
342 
344  virtual size_t getLocalNumCols() const = 0;
345 
347  virtual void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView<LocalOrdinal> &Indices, const ArrayView<Scalar> &Values, size_t &NumEntries) const = 0;
348 
350  virtual bool hasMatrix() const = 0;
351 
353  virtual LocalOrdinal GetStorageBlockSize() const = 0;
354 
359 
360 }; // CrsMatrix class
361 
362 } // namespace Xpetra
363 
364 #define XPETRA_CRSMATRIX_SHORT
365 #endif // XPETRA_CRSMATRIX_HPP
virtual ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const =0
Returns the Frobenius norm of the matrix.
virtual void expertStaticFillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &importer=Teuchos::null, const RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > &exporter=Teuchos::null, const RCP< ParameterList > &params=Teuchos::null)=0
Expert static fill complete.
virtual RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const =0
Returns the CrsGraph associated with this matrix.
virtual global_size_t getGlobalNumRows() const =0
Number of global elements in the row map of this matrix.
virtual void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const =0
Computes the sparse matrix-multivector multiplication.
virtual void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)=0
Replace the diagonal entries of the matrix.
virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const =0
Returns the Map that describes the row distribution in this matrix.
virtual void residual(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const =0
Compute a residual R = B - (*this) * X.
virtual void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Replace matrix entries, using local IDs.
virtual void removeEmptyProcessesInPlace(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)=0
virtual size_t getLocalNumCols() const =0
Returns the number of matrix columns owned on the calling node.
virtual void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Right scale matrix using the given vector entries.
GlobalOrdinal global_ordinal_type
virtual void scale(const Scalar &alpha)=0
Scale the current values of a matrix, this = alpha*this.
virtual bool hasMatrix() const =0
Does this have an underlying matrix.
virtual void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null)=0
Signal that data entry is complete, specifying domain and range maps.
virtual void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Insert matrix entries, using global IDs.
virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const =0
Returns the current number of entries in the specified global row.
virtual size_t getLocalNumEntries() const =0
Returns the local number of entries in this matrix.
virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const =0
Returns the current number of entries on this node in the specified local row.
virtual bool isLocallyIndexed() const =0
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
virtual void setAllToScalar(const Scalar &alpha)=0
Set all matrix entries equal to scalarThis.
LocalOrdinal local_ordinal_type
virtual void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)=0
Sets the 1D pointer arrays of the graph.
virtual void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const =0
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage...
virtual std::string description() const =0
A simple one-line description of this object.
virtual void allocateAllValues(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind, ArrayRCP< Scalar > &values)=0
Allocates and returns ArrayRCPs of the Crs arrays — This is an Xpetra-only routine.
virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const =0
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
virtual void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const =0
Extract a const, non-persisting view of local indices in a specified row of the matrix.
virtual global_size_t getGlobalNumCols() const =0
Number of global columns in the matrix.
Kokkos::StaticCrsGraph< int, Kokkos::LayoutLeft, execution_space, void, size_t > local_graph_type
virtual bool supportsRowViews() const =0
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
virtual void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Replace matrix entries, using global IDs.
virtual void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const =0
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const =0
Print the object with some verbosity level to an FancyOStream object.
virtual void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newDomainMap, Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &newImporter)=0
Replaces the current domainMap and importer with the user-specified objects.
virtual void resumeFill(const RCP< ParameterList > &params=null)=0
virtual size_t getLocalMaxNumRowEntries() const =0
Returns the maximum number of entries across all rows/columns on this node.
virtual local_matrix_type::HostMirror getLocalMatrixHost() const =0
virtual void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Insert matrix entries, using local IDs.
virtual bool haveGlobalConstants() const =0
Returns true if globalConstants have been computed; false otherwise.
typename Kokkos::ArithTraits< double >::val_type impl_scalar_type
size_t global_size_t
Global size_t object.
virtual size_t getLocalNumRows() const =0
Returns the number of matrix rows owned on the calling node.
virtual global_size_t getGlobalNumEntries() const =0
Returns the global number of entries in this matrix.
virtual size_t getGlobalMaxNumRowEntries() const =0
Returns the maximum number of entries across all rows/columns on all nodes.
virtual void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const =0
Get a copy of the diagonal entries owned by this node, with local row indices.
virtual ~CrsMatrix()
Destructor.
virtual void setObjectLabel(const std::string &objectLabel)=0
virtual bool isFillComplete() const =0
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const =0
Returns the Map that describes the column distribution in this matrix.
virtual local_matrix_type getLocalMatrixDevice() const =0
virtual void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const =0
Get offsets of the diagonal entries in the matrix.
virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const =0
Returns the Map associated with the range of this operator, which must be compatible with Y...
virtual bool isGloballyIndexed() const =0
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
KokkosSparse::CrsMatrix< impl_scalar_type, int, execution_space, void, typename local_graph_type::size_type > local_matrix_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
virtual bool isFillActive() const =0
Returns true if the matrix is in edit mode.
virtual LocalOrdinal GetStorageBlockSize() const =0
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
virtual void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const =0
Gets the 1D pointer arrays of the graph.
virtual void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Left scale matrix using the given vector entries.
virtual void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const =0
Extract a const, non-persisting view of global indices in a specified row of the matrix.