Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_CrsMatrixWrap_def.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_CRSMATRIXWRAP_DEF_HPP
13#define XPETRA_CRSMATRIXWRAP_DEF_HPP
14
15#include <Tpetra_KokkosCompat_DefaultNode.hpp>
16
18#include <Teuchos_Hashtable.hpp>
19
20#include "Xpetra_ConfigDefs.hpp"
21#include "Xpetra_Exceptions.hpp"
22
23#include "Xpetra_MultiVector.hpp"
24#include "Xpetra_CrsGraph.hpp"
25#include "Xpetra_CrsMatrix.hpp"
27
28#include "Xpetra_Matrix.hpp"
29
31
32namespace Xpetra {
33template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
39
40template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
42 size_t maxNumEntriesPerRow)
43 : finalDefaultView_(false) {
44 matrixData_ = CrsMatrixFactory::Build(rowMap, maxNumEntriesPerRow);
46}
47
48template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
50 const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc)
51 : finalDefaultView_(false) {
52 matrixData_ = CrsMatrixFactory::Build(rowMap, NumEntriesPerRowToAlloc);
54}
55
56template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
58 : finalDefaultView_(false) {
59 // Set matrix data
60 matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, maxNumEntriesPerRow);
61
62 // Default view
64}
65
66template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
68 : finalDefaultView_(false) {
69 // Set matrix data
70 matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, NumEntriesPerRowToAlloc);
71
72 // Default view
74}
75
76template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
78 : finalDefaultView_(false) {
79 // Set matrix data
80 matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, lclMatrix, params);
81
82 // Default view
84}
85
86template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
88 const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap,
90 : finalDefaultView_(false) {
91 // Set matrix data
92 matrixData_ = CrsMatrixFactory::Build(lclMatrix, rowMap, colMap, domainMap, rangeMap, params);
93
94 // Default view
96}
97
98template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
100 : finalDefaultView_(matrix->isFillComplete()) {
101 // Set matrix data
102 matrixData_ = matrix;
103
104 // Default view
106}
107
108template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
110 : finalDefaultView_(false) {
111 // Set matrix data
112 matrixData_ = CrsMatrixFactory::Build(graph, paramList);
113
114 // Default view
116}
117
118template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
121
122 const RCP<ParameterList> &paramList)
123 : finalDefaultView_(false) {
124 // Set matrix data
125 matrixData_ = CrsMatrixFactory::Build(graph, values, paramList);
126
127 // Default view
129}
130
131template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
133
134template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
136 matrixData_->insertGlobalValues(globalRow, cols, vals);
137}
138
139template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
141 matrixData_->insertLocalValues(localRow, cols, vals);
142}
143
144template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
147 const ArrayView<const Scalar> &vals) { matrixData_->replaceGlobalValues(globalRow, cols, vals); }
148
149template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
152 const ArrayView<const Scalar> &vals) { matrixData_->replaceLocalValues(localRow, cols, vals); }
153
154template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
155void CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::setAllToScalar(const Scalar &alpha) { matrixData_->setAllToScalar(alpha); }
156
157template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
159 matrixData_->scale(alpha);
160}
161
162template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
164 matrixData_->resumeFill(params);
165}
166
167template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
169 matrixData_->fillComplete(domainMap, rangeMap, params);
170
171 // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
172 updateDefaultView();
173}
174
175template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
177 matrixData_->fillComplete(params);
178
179 // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
180 updateDefaultView();
181}
182
183template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
185 return matrixData_->getGlobalNumRows();
186}
187
188template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
190 return matrixData_->getGlobalNumCols();
191}
192
193template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
195 return matrixData_->getLocalNumRows();
196}
197
198template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
200 return matrixData_->getGlobalNumEntries();
201}
202
203template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
205 return matrixData_->getLocalNumEntries();
206}
207
208template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
210 return matrixData_->getNumEntriesInLocalRow(localRow);
211}
212
213template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
215 return matrixData_->getNumEntriesInGlobalRow(globalRow);
216}
217
218template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
220 return matrixData_->getGlobalMaxNumRowEntries();
221}
222
223template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
225 return matrixData_->getLocalMaxNumRowEntries();
226}
227
228template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
230 return matrixData_->isLocallyIndexed();
231}
232
233template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
235 return matrixData_->isGloballyIndexed();
236}
237
238template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
240 return matrixData_->isFillComplete();
241}
242
243template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
245 const ArrayView<LocalOrdinal> &Indices,
246 const ArrayView<Scalar> &Values,
247 size_t &NumEntries) const {
248 matrixData_->getLocalRowCopy(LocalRow, Indices, Values, NumEntries);
249}
250
251template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
253 matrixData_->getGlobalRowView(GlobalRow, indices, values);
254}
255
256template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
258 matrixData_->getLocalRowView(LocalRow, indices, values);
259}
260
261template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
265
266template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
268 matrixData_->getLocalDiagOffsets(offsets);
269}
270
271template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
275
276template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
280
281template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
285
286template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
290
291template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
293 return matrixData_->haveGlobalConstants();
294}
295
296template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
304
305template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
308 Teuchos::ETransp mode,
309 Scalar alpha,
310 Scalar beta,
311 bool sumInterfaceValues,
312 const RCP<Import<LocalOrdinal, GlobalOrdinal, Node>> &regionInterfaceImporter,
313 const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const {
314 matrixData_->apply(X, Y, mode, alpha, beta, sumInterfaceValues, regionInterfaceImporter, regionInterfaceLIDs);
315}
316
317template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
321
322template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
326
327template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
329
330template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
332 TEUCHOS_TEST_FOR_EXCEPTION(Matrix::operatorViewTable_.containsKey(viewLabel) == false, Xpetra::Exceptions::RuntimeError, "Xpetra::Matrix.GetColMap(): view '" + viewLabel + "' does not exist.");
333 updateDefaultView(); // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated.
334 return Matrix::operatorViewTable_.get(viewLabel)->GetColMap();
335}
336
337template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
339 matrixData_->removeEmptyProcessesInPlace(newMap);
340 this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetRowMap(matrixData_->getRowMap());
341 this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetColMap(matrixData_->getColMap());
342}
343
344template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
348
349template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
352 const CrsMatrixWrap &sourceWrp = dynamic_cast<const CrsMatrixWrap &>(source);
353 matrixData_->doImport(*sourceWrp.getCrsMatrix(), importer, CM);
354}
355
356template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
359 const CrsMatrixWrap &destWrp = dynamic_cast<const CrsMatrixWrap &>(dest);
360 matrixData_->doExport(*destWrp.getCrsMatrix(), importer, CM);
361}
362
363template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
366 const CrsMatrixWrap &sourceWrp = dynamic_cast<const CrsMatrixWrap &>(source);
367 matrixData_->doImport(*sourceWrp.getCrsMatrix(), exporter, CM);
368}
369
370template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
373 const CrsMatrixWrap &destWrp = dynamic_cast<const CrsMatrixWrap &>(dest);
374 matrixData_->doExport(*destWrp.getCrsMatrix(), exporter, CM);
375}
376
377template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
379 return "Xpetra::CrsMatrixWrap";
380}
381
382template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
384 // Teuchos::EVerbosityLevel vl = verbLevel;
385 // if (vl == VERB_DEFAULT) vl = VERB_LOW;
386 // RCP<const Comm<int> > comm = this->getComm();
387 // const int myImageID = comm->getRank(),
388 // numImages = comm->getSize();
389
390 // if (myImageID == 0) out << this->description() << std::endl;
391
392 matrixData_->describe(out, verbLevel);
393
394 // Teuchos::OSTab tab(out);
395}
396
397template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
400 matrixData_->setObjectLabel(objectLabel);
401}
402
403template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
408template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
411 return matrixData_->getLocalMatrixDevice();
412}
413
414template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
416
417template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
419
420template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
422
423// Default view is created after fillComplete()
424// Because ColMap might not be available before fillComplete().
425template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
427 // Create default view
428 this->defaultViewLabel_ = "point";
429 this->CreateView(this->GetDefaultViewLabel(), matrixData_->getRowMap(), matrixData_->getColMap());
430
431 // Set current view
432 this->currentViewLabel_ = this->GetDefaultViewLabel();
433}
434
435template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
437 if ((finalDefaultView_ == false) && matrixData_->isFillComplete()) {
438 // Update default view with the colMap
439 Matrix::operatorViewTable_.get(Matrix::GetDefaultViewLabel())->SetColMap(matrixData_->getColMap());
440 finalDefaultView_ = true;
441 }
442}
443
444// Expert only
445template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
447 // Clear the old view table
449 Matrix::operatorViewTable_ = dummy_table;
450
451 finalDefaultView_ = M->isFillComplete();
452 // Set matrix data
453 matrixData_ = M;
454
455 // Default view
456 CreateDefaultView();
457}
458
459template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
461 return matrixData_->GetStorageBlockSize();
462}
463
464template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
471
472} // namespace Xpetra
473
474#endif // ifndef XPETRA_CRSMATRIXWRAP_DEF_HPP
virtual void setObjectLabel(const std::string &objectLabel)
static RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap)
Constructor for empty matrix (intended use is an import/export target - can't insert entries directly...
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...
Concrete implementation of Xpetra::Matrix.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale matrix using the given vector entries.
LocalOrdinal GetStorageBlockSize() const
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
void residual(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const
Compute a residual R = B - (*this) * X.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the specified global row.
RCP< const CrsGraph > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< Teuchos::ParameterList > &params=null)
Signal that data entry is complete, specifying domain and range maps.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale matrix using the given vector entries.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Get Frobenius norm of the matrix.
std::string description() const
Return a simple one-line description of this object.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
bool hasCrsGraph() const
Supports the getCrsGraph() call.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type local_matrix_type
void resumeFill(const RCP< ParameterList > &params=null)
global_size_t getGlobalNumCols() const
Returns the number of global columns in the matrix.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
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.
global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
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 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.
void setObjectLabel(const std::string &objectLabel)
virtual local_matrix_type::host_mirror_type getLocalMatrixHost() const
const RCP< const Map > & getColMap() const
Returns the Map that describes the column distribution in this matrix. This might be null until fillC...
const RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
virtual local_matrix_type getLocalMatrixDevice() const
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void replaceCrsMatrix(RCP< CrsMatrix > &M)
Expert only.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void doExport(const Matrix &dest, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void getLocalDiagCopy(Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row idices.
virtual ~CrsMatrixWrap()
Destructor.
const RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
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.
size_t getLocalNumRows() const
Returns the number of matrix rows owned on the calling node.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
CrsMatrixWrap(const RCP< const Map > &rowMap)
Constructor for a dynamic profile matrix (Epetra only)
void doImport(const Matrix &source, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
RCP< CrsMatrix > getCrsMatrix() const
virtual void apply(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode.
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
Exception throws to report errors in the internal logical of the program.
Xpetra-specific matrix class.
const viewLabel_t & GetCurrentViewLabel() const
const viewLabel_t & GetDefaultViewLabel() const
Teuchos::Hashtable< viewLabel_t, RCP< MatrixView > > operatorViewTable_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
size_t global_size_t
Global size_t object.
CombineMode
Xpetra::Combine Mode enumerable type.