Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_BlockedCrsMatrix_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#ifndef XPETRA_BLOCKEDCRSMATRIX_DECL_HPP
11#define XPETRA_BLOCKEDCRSMATRIX_DECL_HPP
12
13#include <Tpetra_KokkosCompat_DefaultNode.hpp>
14
16#include <Teuchos_Hashtable.hpp>
17
18#include "Xpetra_ConfigDefs.hpp"
19#include "Xpetra_Exceptions.hpp"
20
21#include "Xpetra_MapFactory.hpp"
22#include "Xpetra_MultiVector.hpp"
23#include "Xpetra_BlockedMultiVector.hpp"
24#include "Xpetra_MultiVectorFactory.hpp"
25#include "Xpetra_BlockedVector.hpp"
26#include "Xpetra_CrsGraph.hpp"
27#include "Xpetra_CrsMatrix.hpp"
29
30#include "Xpetra_MapExtractor.hpp"
32
33#include "Xpetra_Matrix.hpp"
34#include "Xpetra_MatrixFactory.hpp"
35#include "Xpetra_CrsMatrixWrap.hpp"
36
37#ifdef HAVE_XPETRA_THYRA
38#include <Thyra_ProductVectorSpaceBase.hpp>
39#include <Thyra_VectorSpaceBase.hpp>
40#include <Thyra_LinearOpBase.hpp>
41#include <Thyra_BlockedLinearOpBase.hpp>
42#include <Thyra_PhysicallyBlockedLinearOpBase.hpp>
43#include "Xpetra_ThyraUtils.hpp"
44#endif
45
46#include "Xpetra_VectorFactory.hpp"
47
52namespace Xpetra {
53
54#ifdef HAVE_XPETRA_THYRA
55template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
56class ThyraUtils;
57#endif
58
59typedef std::string viewLabel_t;
60
61template <class Scalar,
62 class LocalOrdinal,
63 class GlobalOrdinal,
64 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
65class BlockedCrsMatrix : public Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
66 public:
67 typedef Scalar scalar_type;
68 typedef LocalOrdinal local_ordinal_type;
69 typedef GlobalOrdinal global_ordinal_type;
70 typedef Node node_type;
71
72 private:
73#undef XPETRA_BLOCKEDCRSMATRIX_SHORT
75
76 public:
78
79
81
87 const Teuchos::RCP<const BlockedMap>& domainMaps,
88 size_t numEntriesPerRow);
89
91
99 Teuchos::RCP<const MapExtractor>& domainMapExtractor,
100 size_t numEntriesPerRow);
101
102#ifdef HAVE_XPETRA_THYRA
104
109 BlockedCrsMatrix(const Teuchos::RCP<const Thyra::BlockedLinearOpBase<Scalar> >& thyraOp, const Teuchos::RCP<const Teuchos::Comm<int> >& /* comm */);
110
111 private:
113
121
122 public:
123#endif
124
127
129
131
132
134
159 void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal>& cols, const ArrayView<const Scalar>& vals);
160
162
172 void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal>& cols, const ArrayView<const Scalar>& vals);
173
175
177
185 void replaceGlobalValues(GlobalOrdinal globalRow,
187 const ArrayView<const Scalar>& vals);
188
190
194 void replaceLocalValues(LocalOrdinal localRow,
196 const ArrayView<const Scalar>& vals);
197
199 virtual void setAllToScalar(const Scalar& alpha);
200
202 void scale(const Scalar& alpha);
203
205
207
208
216 void resumeFill(const RCP<ParameterList>& params = null);
217
223 void fillComplete(const RCP<const Map>& domainMap, const RCP<const Map>& rangeMap, const RCP<ParameterList>& params = null);
224
239 void fillComplete(const RCP<ParameterList>& params = null);
240
242
244
247
249
252
254 size_t getLocalNumRows() const;
255
258
260 size_t getLocalNumEntries() const;
261
263
264 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
265
267
268 size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
269
271
273 size_t getGlobalMaxNumRowEntries() const;
274
276
278 size_t getLocalMaxNumRowEntries() const;
279
281
284 bool isLocallyIndexed() const;
285
287
290 bool isGloballyIndexed() const;
291
293 bool isFillComplete() const;
294
296
310 virtual void getLocalRowCopy(LocalOrdinal LocalRow,
311 const ArrayView<LocalOrdinal>& Indices,
312 const ArrayView<Scalar>& Values,
313 size_t& NumEntries) const;
314
316
325 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal>& indices, ArrayView<const Scalar>& values) const;
326
328
337 void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal>& indices, ArrayView<const Scalar>& values) const;
338
340
343 void getLocalDiagCopy(Vector& diag) const;
344
346 void leftScale(const Vector& x);
347
349 void rightScale(const Vector& x);
350
353
355 virtual bool haveGlobalConstants() const;
356
358
360
361
363
383
385
386
388 virtual void apply(const MultiVector& X, MultiVector& Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues,
389 const RCP<Xpetra::Import<LocalOrdinal, GlobalOrdinal, Node> >& regionInterfaceImporter,
390 const Teuchos::ArrayRCP<LocalOrdinal>& regionInterfaceLIDs) const;
391
393
396 virtual void apply(const MultiVector& X, MultiVector& Y,
398 Scalar alpha = ScalarTraits<Scalar>::one(),
399 Scalar beta = ScalarTraits<Scalar>::zero()) const;
400
403
406
408 const RCP<const Map> getDomainMap() const;
409
411 RCP<const Map> getDomainMap(size_t i) const;
412
414 RCP<const Map> getDomainMap(size_t i, bool bThyraMode) const;
415
418
421
423 const RCP<const Map> getRangeMap() const;
424
426 RCP<const Map> getRangeMap(size_t i) const;
427
429 RCP<const Map> getRangeMap(size_t i, bool bThyraMode) const;
430
433
436
438
440 //{@
441
450 virtual void bgs_apply(
451 const MultiVector& X,
452 MultiVector& Y,
453 size_t row,
455 Scalar alpha = ScalarTraits<Scalar>::one(),
456 Scalar beta = ScalarTraits<Scalar>::zero()
457 ) const;
458
460
462 //{@
463
465 const Teuchos::RCP<const Map> getMap() const;
466
468 void doImport(const Matrix& source, const Import& importer, CombineMode CM);
469
471 void doExport(const Matrix& dest, const Import& importer, CombineMode CM);
472
474 void doImport(const Matrix& source, const Export& exporter, CombineMode CM);
475
477 void doExport(const Matrix& dest, const Export& exporter, CombineMode CM);
478
479 // @}
480
482
483
485 std::string description() const;
486
489
491
492 void setObjectLabel(const std::string& objectLabel);
494
496 bool hasCrsGraph() const;
497
500
502
504
505
506 virtual bool isDiagonal() const;
507
509 virtual size_t Rows() const;
510
512 virtual size_t Cols() const;
513
516
519
521 Teuchos::RCP<Matrix> getMatrix(size_t r, size_t c) const;
522
524 // void setMatrix(size_t r, size_t c, Teuchos::RCP<CrsMatrix> mat) {
525 void setMatrix(size_t r, size_t c, Teuchos::RCP<Matrix> mat);
526
528 // NOTE: This is a rather expensive operation, since all blocks are copied
529 // into a new big CrsMatrix
532
537#if KOKKOS_VERSION >= 40799
538 typename local_matrix_type::host_mirror_type getLocalMatrixHost() const;
539#else
540 typename local_matrix_type::HostMirror getLocalMatrixHost() const;
541#endif
542
543#ifdef HAVE_XPETRA_THYRA
545#endif
547 LocalOrdinal GetStorageBlockSize() const;
548
550 void residual(const MultiVector& X,
551 const MultiVector& B,
552 MultiVector& R) const;
553
554 private:
557
559
569 void Add(const Matrix& A, const Scalar scalarA, Matrix& B, const Scalar scalarB) const;
570
572
573 // Default view is created after fillComplete()
574 // Because ColMap might not be available before fillComplete().
575 void CreateDefaultView();
576
577 private:
581
582 std::vector<Teuchos::RCP<Matrix> > blocks_;
583#ifdef HAVE_XPETRA_THYRA
585#endif
588};
589
590} // namespace Xpetra
591
592#define XPETRA_BLOCKEDCRSMATRIX_SHORT
593#endif /* XPETRA_BLOCKEDCRSMATRIX_DECL_HPP */
static const EVerbosityLevel verbLevel_default
Teuchos::RCP< const MapExtractor > domainmaps_
full domain map together with all partial domain maps
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
virtual ~BlockedCrsMatrix()
Destructor.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
virtual size_t Cols() const
number of column blocks
bool isLocallyIndexed() const
If matrix indices of all matrix blocks are in the local range, this function returns true....
Teuchos::RCP< const MapExtractor > rangemaps_
full range map together with all partial domain maps
const RCP< const Map > getRangeMap() const
Returns the Map associated with the range of this operator.
void setMatrix(size_t r, size_t c, Teuchos::RCP< Matrix > mat)
set matrix block
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
size_t getLocalNumRows() const
Returns the number of matrix rows owned on the calling node.
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null)
Signal that data entry is complete.
bool bDomainThyraMode_
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub...
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
const Teuchos::RCP< const Map > getMap() const
Implements DistObject interface.
bool hasCrsGraph() const
Supports the getCrsGraph() call.
Teuchos::RCP< Matrix > getMatrix(size_t r, size_t c) const
return block (r,c)
std::vector< Teuchos::RCP< Matrix > > blocks_
row major matrix block storage
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
RCP< const BlockedMap > getBlockedDomainMap() const
Returns the BlockedMap associated with the domain of this operator.
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode.
Teuchos::RCP< Matrix > Merge() const
merge BlockedCrsMatrix blocks in a CrsMatrix
local_matrix_type::HostMirror getLocalMatrixHost() const
Access the underlying local Kokkos::CrsMatrix object.
void leftScale(const Vector &x)
Left scale matrix using the given vector entries.
const RCP< const Map > getDomainMap() const
Returns the Map associated with the domain of this operator.
Teuchos::RCP< Matrix > getCrsMatrix() const
return unwrap 1x1 blocked operators
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
void doImport(const Matrix &source, const Import &importer, CombineMode CM)
Import.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the matrix.
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
bool is_diagonal_
If we're diagonal, a bunch of the extraction stuff should work.
RCP< const CrsGraph > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
virtual bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
void Add(const Matrix &A, const Scalar scalarA, Matrix &B, const Scalar scalarB) const
Add a Xpetra::CrsMatrix to another: B = B*scalarB + A*scalarA.
bool bRangeThyraMode_
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub...
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
virtual size_t Rows() const
number of row blocks
RCP< const MapExtractor > getRangeMapExtractor() const
Returns map extractor class for range map.
RCP< const BlockedMap > getBlockedRangeMap() const
Returns the BlockedMap associated with the range of this operator.
LocalOrdinal GetStorageBlockSize() const
Returns the block size of the storage mechanism.
virtual void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP< Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > > &regionInterfaceImporter, const Teuchos::ArrayRCP< LocalOrdinal > &regionInterfaceLIDs) const
sparse matrix-multivector multiplication for the region layout matrices (currently no blocked impleme...
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the specified (locally owned) global row.
CrsMatrix::local_matrix_type local_matrix_type
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
virtual ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Get Frobenius norm of the matrix.
virtual void bgs_apply(const MultiVector &X, MultiVector &Y, size_t row, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Special multiplication routine (for BGS/Jacobi smoother)
RCP< const Map > getFullDomainMap() const
Returns the Map associated with the full domain of this operator.
void rightScale(const Vector &x)
Right scale matrix using the given vector entries.
void setObjectLabel(const std::string &objectLabel)
RCP< const Map > getFullRangeMap() const
Returns the Map associated with the full range of this operator.
void resumeFill(const RCP< ParameterList > &params=null)
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
void getLocalDiagCopy(Vector &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
Teuchos::RCP< Matrix > getInnermostCrsMatrix()
helper routine recursively returns the first inner-most non-null matrix block from a (nested) blocked...
global_size_t getGlobalNumRows() const
Returns the number of global rows.
RCP< const MapExtractor > getDomainMapExtractor() const
Returns map extractor for domain map.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
virtual void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
void doExport(const Matrix &dest, const Import &importer, CombineMode CM)
Export.
std::string description() const
Return a simple one-line description of this object.
local_matrix_type getLocalMatrixDevice() const
Access the underlying local Kokkos::CrsMatrix object.
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...
Xpetra-specific matrix class.
size_t global_size_t
Global size_t object.
CombineMode
Xpetra::Combine Mode enumerable type.