42#ifndef TPETRA_REPLACEDIAGONALCRSMATRIX_DEF_HPP
43#define TPETRA_REPLACEDIAGONALCRSMATRIX_DEF_HPP
49#include "Tpetra_CrsMatrix.hpp"
50#include "Tpetra_Vector.hpp"
54template<
class SC,
class LO,
class GO,
class NT>
63 const LO oneLO = Teuchos::OrdinalTraits<LO>::one();
66 LO numReplacedDiagEntries = 0;
70 if (rowMapPtr.is_null() || rowMapPtr->getComm().is_null()) {
74 return numReplacedDiagEntries;
78 TEUCHOS_TEST_FOR_EXCEPTION
79 (colMapPtr.get () ==
nullptr, std::invalid_argument,
80 "replaceDiagonalCrsMatrix: "
81 "Input matrix must have a nonnull column Map.");
89 TEUCHOS_TEST_FOR_EXCEPTION(!rowMap.isSameAs(*newDiag.
getMap()), std::runtime_error,
90 "Row map of matrix and map of input vector do not match.");
96 typename crs_matrix_type::execution_space().fence();
98 if (isFillCompleteOnInput)
101 Teuchos::ArrayRCP<const SC> newDiagData = newDiag.
getVector(0)->getData();
102 LO numReplacedEntriesPerRow = 0;
104 auto invalid = Teuchos::OrdinalTraits<LO>::invalid();
107 for (LO lclRowInd = 0; lclRowInd < myNumRows; ++lclRowInd) {
110 const GO gblInd = rowMap.getGlobalElement(lclRowInd);
111 const LO lclColInd = colMap.getLocalElement(gblInd);
115 if (lclColInd == invalid)
continue;
117 const SC vals[] = {
static_cast<SC
>(newDiagData[lclRowInd])};
118 const LO cols[] = {lclColInd};
130 if (numReplacedEntriesPerRow == oneLO) {
131 ++numReplacedDiagEntries;
135 if (isFillCompleteOnInput)
138 return numReplacedDiagEntries;
148#define TPETRA_REPLACEDIAGONALCRSMATRIX_INSTANT(SCALAR,LO,GO,NODE) \
151 LO replaceDiagonalCrsMatrix( \
152 CrsMatrix<SCALAR, LO, GO, NODE>& matrix, \
153 const Vector<SCALAR, LO, GO, NODE>& newDiag); \
Declaration of Tpetra::Details::Behavior, a class that describes Tpetra's behavior.
Sparse matrix that presents a row-oriented interface that lets users read or modify entries.
Teuchos::RCP< const map_type > getColMap() const override
The Map that describes the column distribution in this matrix.
void fillComplete(const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Tell the matrix that you are done changing its structure or values, and that you are ready to do comp...
size_t getLocalNumRows() const override
The number of matrix rows owned by the calling process.
bool isFillComplete() const override
Whether the matrix is fill complete.
Teuchos::RCP< const map_type > getRowMap() const override
The Map that describes the row distribution in this matrix.
local_ordinal_type replaceLocalValues(const local_ordinal_type localRow, const Kokkos::View< const local_ordinal_type *, Kokkos::AnonymousSpace > &inputInds, const Kokkos::View< const impl_scalar_type *, Kokkos::AnonymousSpace > &inputVals)
Replace one or more entries' values, using local row and column indices.
void resumeFill(const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Resume operations that may change the values or structure of the matrix.
static bool debug()
Whether Tpetra is in debug mode.
virtual Teuchos::RCP< const map_type > getMap() const
The Map describing the parallel distribution of this object.
A parallel distribution of indices over processes.
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(const size_t j) const
Return a Vector which is a const view of column j.
A distributed dense vector.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
LO replaceDiagonalCrsMatrix(::Tpetra::CrsMatrix< SC, LO, GO, NT > &matrix, const ::Tpetra::Vector< SC, LO, GO, NT > &newDiag)
Replace diagonal entries of an input Tpetra::CrsMatrix matrix with values given in newDiag.