Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_Matrix_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Xpetra: A linear algebra interface package
4//
5// Copyright 2012 NTESS and the Xpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10// WARNING: This code is experimental. Backwards compatibility should not be expected.
11
12#ifndef XPETRA_MATRIX_DECL_HPP
13#define XPETRA_MATRIX_DECL_HPP
14
15#include <Tpetra_KokkosCompat_DefaultNode.hpp>
16
17#include "Xpetra_ConfigDefs.hpp"
18#include "Xpetra_Exceptions.hpp"
19
20#include "Xpetra_MultiVector.hpp"
21#include "Xpetra_CrsGraph.hpp"
22#include "Xpetra_CrsMatrix.hpp"
24#include "Xpetra_MatrixView.hpp"
25#include "Xpetra_Operator.hpp"
26#include "Xpetra_StridedMap.hpp"
27#include "Xpetra_StridedMapFactory.hpp"
28
30#include <Teuchos_Hashtable.hpp>
31
36namespace Xpetra {
37
54typedef std::string viewLabel_t;
55
56template <class Scalar,
57 class LocalOrdinal,
58 class GlobalOrdinal,
59 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
60class Matrix : public Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
66
67 public:
68 typedef Scalar scalar_type;
69 typedef LocalOrdinal local_ordinal_type;
70 typedef GlobalOrdinal global_ordinal_type;
71 typedef Node node_type;
72
73#ifdef HAVE_XPETRA_TPETRA
77#endif
78
79#if KOKKOS_VERSION >= 40799
80 using ATS = KokkosKernels::ArithTraits<Scalar>;
81#else
82 using ATS = Kokkos::ArithTraits<Scalar>;
83#endif
84 using impl_scalar_type = typename ATS::val_type;
85
87
88
89 Matrix();
90
92 virtual ~Matrix();
93
95
97
98 void CreateView(viewLabel_t viewLabel, const RCP<const Map> &rowMap, const RCP<const Map> &colMap);
99
100 // JG TODO: why this is a member function??
101 void CreateView(const viewLabel_t viewLabel, const RCP<const Matrix> &A, bool transposeA = false, const RCP<const Matrix> &B = Teuchos::null, bool transposeB = false);
102
104 void PrintViews(Teuchos::FancyOStream &out) const;
105
106 void RemoveView(const viewLabel_t viewLabel);
107
108 const viewLabel_t SwitchToView(const viewLabel_t viewLabel);
109
110 bool IsView(const viewLabel_t viewLabel) const;
111
113
114 const viewLabel_t &GetDefaultViewLabel() const;
115
116 const viewLabel_t &GetCurrentViewLabel() const;
117
119
121
122
124
135 virtual void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals) = 0;
136
138
145 virtual void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals) = 0;
146
148
153 virtual void replaceGlobalValues(GlobalOrdinal globalRow,
155 const ArrayView<const Scalar> &vals) = 0;
156
158
161 virtual void replaceLocalValues(LocalOrdinal localRow,
163 const ArrayView<const Scalar> &vals) = 0;
164
166 virtual void setAllToScalar(const Scalar &alpha) = 0;
167
169 virtual void scale(const Scalar &alpha) = 0;
170
172
174
175
184 virtual void resumeFill(const RCP<ParameterList> &params = null) = 0;
185
197 virtual void fillComplete(const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap, const RCP<ParameterList> &params = null) = 0;
198
212 // TODO : Get ride of "Tpetra"::OptimizeOption
213 virtual void fillComplete(const RCP<ParameterList> &params = null) = 0;
214
216
218
219
221 virtual const RCP<const Map> &getRowMap() const;
222
224 virtual const RCP<const Map> &getRowMap(viewLabel_t viewLabel) const;
225
228 virtual const RCP<const Map> &getColMap() const;
229
231 virtual const RCP<const Map> &getColMap(viewLabel_t viewLabel) const;
232
234
236 virtual global_size_t getGlobalNumRows() const = 0;
237
239
241 virtual global_size_t getGlobalNumCols() const = 0;
242
244 virtual size_t getLocalNumRows() const = 0;
245
248
250 virtual size_t getLocalNumEntries() const = 0;
251
253
254 virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const = 0;
255
257
258 virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const = 0;
259
261
263 virtual size_t getGlobalMaxNumRowEntries() const = 0;
264
266
268 virtual size_t getLocalMaxNumRowEntries() const = 0;
269
271 virtual bool isLocallyIndexed() const = 0;
272
274 virtual bool isGloballyIndexed() const = 0;
275
277 virtual bool isFillComplete() const = 0;
278
280
292 virtual void getLocalRowCopy(LocalOrdinal LocalRow,
293 const ArrayView<LocalOrdinal> &Indices,
294 const ArrayView<Scalar> &Values,
295 size_t &NumEntries) const = 0;
296
298
307 virtual void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) const = 0;
308
310
319 virtual void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) const = 0;
320
322
325
328
331
334
336 virtual bool haveGlobalConstants() const = 0;
337
339
341 //{@
342
345
346 // TODO: first argument of doImport/doExport should be a Xpetra::DistObject
347
349 virtual void doImport(const Matrix &source,
351
353 virtual void doExport(const Matrix &dest,
355
357 virtual void doImport(const Matrix &source,
359
361 virtual void doExport(const Matrix &dest,
363
364 // @}
365
367
368
370
371
373 virtual std::string description() const = 0;
374
378
380
381 virtual void setObjectLabel(const std::string &objectLabel) = 0;
383
384 // JG: Added:
385
387 virtual bool hasCrsGraph() const = 0;
388
390 virtual RCP<const CrsGraph> getCrsGraph() const = 0;
391
392 // To keep the normal virtual matrix-multivector definition of apply before overloading with the region variant
393 using Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>::apply;
394
398 Teuchos::ETransp mode,
399 Scalar alpha,
400 Scalar beta,
401 bool sumInterfaceValues,
402 const RCP<Import<LocalOrdinal, GlobalOrdinal, Node> > &regionInterfaceImporter,
403 const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const = 0;
404
405 // ----------------------------------------------------------------------------------
406 // "TEMPORARY" VIEW MECHANISM
413 void SetFixedBlockSize(LocalOrdinal blksize, GlobalOrdinal offset = 0);
414
415 //==========================================================================
416
417 LocalOrdinal GetFixedBlockSize() const; // TODO: why LocalOrdinal?
418
420 bool IsFixedBlockSizeSet() const;
421
423 virtual LocalOrdinal GetStorageBlockSize() const = 0;
424
425 // ----------------------------------------------------------------------------------
426
427 virtual void SetMaxEigenvalueEstimate(Scalar const &sigma);
428
429 // ----------------------------------------------------------------------------------
430
431 virtual Scalar GetMaxEigenvalueEstimate() const;
432
433 // ----------------------------------------------------------------------------------
434#ifdef HAVE_XPETRA_TPETRA
436#if KOKKOS_VERSION >= 40799
437 virtual typename local_matrix_type::host_mirror_type getLocalMatrixHost() const = 0;
438#else
439 virtual typename local_matrix_type::HostMirror getLocalMatrixHost() const = 0;
440#endif
441#else
442#ifdef __GNUC__
443#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."
444#endif
445#endif
446 // ----------------------------------------------------------------------------------
447
452
453 protected:
454 Teuchos::Hashtable<viewLabel_t, RCP<MatrixView> > operatorViewTable_; // hashtable storing the operator views (keys = view names, values = views).
455
456 viewLabel_t defaultViewLabel_; // label of the view associated with inital Matrix construction
457 viewLabel_t currentViewLabel_; // label of the current view
458
459}; // class Matrix
460
461} // namespace Xpetra
462
463#define XPETRA_MATRIX_SHORT
464#endif // XPETRA_MATRIX_DECL_HPP
static const EVerbosityLevel verbLevel_default
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, 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...
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, execution_space, void, typename local_graph_type::size_type > local_matrix_device_type
typename local_matrix_device_type::HostMirror local_matrix_host_type
Xpetra-specific matrix class.
virtual global_size_t getGlobalNumEntries() const =0
Returns the global number of entries in this matrix.
virtual bool haveGlobalConstants() const =0
Returns true if globalConstants have been computed; false otherwise.
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
const viewLabel_t & GetCurrentViewLabel() const
void CreateView(viewLabel_t viewLabel, const RCP< const Map > &rowMap, const RCP< const Map > &colMap)
virtual size_t getGlobalMaxNumRowEntries() const =0
Returns the maximum number of entries across all rows/columns on all nodes.
virtual global_size_t getGlobalNumRows() const =0
Returns the number of global rows in this matrix.
virtual bool hasCrsGraph() const =0
Supports the getCrsGraph() call.
viewLabel_t currentViewLabel_
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > CrsGraph
virtual void SetMaxEigenvalueEstimate(Scalar const &sigma)
Xpetra::CrsMatrixFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixFactory
LocalOrdinal local_ordinal_type
void SetFixedBlockSize(LocalOrdinal blksize, GlobalOrdinal offset=0)
virtual void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Replace matrix entries, using local IDs.
virtual ~Matrix()
Destructor.
virtual const RCP< const Map > & getRowMap() const
Returns the Map that describes the row distribution in this matrix.
virtual void resumeFill(const RCP< ParameterList > &params=null)=0
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type local_matrix_type
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 bool isLocallyIndexed() const =0
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
virtual bool isGloballyIndexed() const =0
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
virtual void doImport(const Matrix &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)=0
Import (using an Exporter).
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
virtual bool isFillComplete() const =0
Returns true if fillComplete() has been called and the matrix is in compute mode.
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 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 void setAllToScalar(const Scalar &alpha)=0
Set all matrix entries equal to scalar.
virtual void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Insert matrix entries, using local IDs.
virtual const RCP< const Map > & getColMap() const
Returns the Map that describes the column distribution in this matrix. This might be null until fillC...
const viewLabel_t SwitchToView(const viewLabel_t viewLabel)
virtual std::string description() const =0
Return a simple one-line description of this object.
virtual size_t getLocalNumRows() const =0
Returns the number of matrix rows owned on the calling node.
virtual void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null)=0
Signal that data entry is complete, specifying domain and range maps.
const viewLabel_t & GetDefaultViewLabel() const
virtual void doExport(const Matrix &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)=0
Export (using an Importer).
viewLabel_t defaultViewLabel_
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_host_type local_matrix_host_type
virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const =0
Returns the current number of entries on this node in the specified local row.
virtual local_matrix_type::HostMirror getLocalMatrixHost() const =0
bool IsView(const viewLabel_t viewLabel) const
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.
GlobalOrdinal global_ordinal_type
virtual void fillComplete(const RCP< ParameterList > &params=null)=0
Signal that data entry is complete.
bool IsFixedBlockSizeSet() const
Returns true, if SetFixedBlockSize has been called before.
virtual void setObjectLabel(const std::string &objectLabel)=0
void RemoveView(const viewLabel_t viewLabel)
virtual size_t getLocalNumEntries() const =0
Returns the local number of entries in this matrix.
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.
const viewLabel_t SwitchToDefaultView()
virtual void doImport(const Matrix &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)=0
Import.
virtual void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Left scale matrix using the given vector entries.
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
Computes the matrix-multivector multiplication for region layout matrices.
virtual void doExport(const Matrix &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)=0
Export.
Xpetra::MatrixView< Scalar, LocalOrdinal, GlobalOrdinal, Node > MatrixView
virtual local_matrix_type getLocalMatrixDevice() const =0
virtual void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Replace matrix entries, using global IDs.
void PrintViews(Teuchos::FancyOStream &out) const
Print all of the views associated with the Matrix.
virtual size_t getLocalMaxNumRowEntries() const =0
Returns the maximum number of entries across all rows/columns on this node.
Teuchos::Hashtable< viewLabel_t, RCP< MatrixView > > operatorViewTable_
virtual LocalOrdinal GetStorageBlockSize() const =0
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
virtual RCP< const CrsGraph > getCrsGraph() const =0
Returns the CrsGraph associated with this matrix.
virtual ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const =0
Get Frobenius norm of the matrix.
LocalOrdinal GetFixedBlockSize() const
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 rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Right scale matrix using the given vector entries.
virtual global_size_t getGlobalNumCols() const =0
Returns the number of global columns in the matrix.
typename ATS::val_type impl_scalar_type
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_device_type local_matrix_device_type
virtual void scale(const Scalar &alpha)=0
Scale the current values of a matrix, this = alpha*this.
Kokkos::ArithTraits< Scalar > ATS
virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const =0
Returns the current number of entries in the specified global row.
virtual const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const =0
Implements DistObject interface.
virtual Scalar GetMaxEigenvalueEstimate() const
virtual void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Insert matrix entries, using global IDs.
size_t global_size_t
Global size_t object.
CombineMode
Xpetra::Combine Mode enumerable type.