Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_replaceDiagonalCrsMatrix_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef TPETRA_REPLACEDIAGONALCRSMATRIX_DEF_HPP
11#define TPETRA_REPLACEDIAGONALCRSMATRIX_DEF_HPP
12
15
17#include "Tpetra_CrsMatrix.hpp"
18#include "Tpetra_Vector.hpp"
19
20namespace Tpetra {
21
22template <class SC, class LO, class GO, class NT>
23LO replaceDiagonalCrsMatrix(CrsMatrix<SC, LO, GO, NT>& matrix,
24 const Vector<SC, LO, GO, NT>& newDiag) {
25 using map_type = Map<LO, GO, NT>;
26 using crs_matrix_type = CrsMatrix<SC, LO, GO, NT>;
27 // using vec_type = ::Tpetra::Vector<SC, LO, GO, NT>;
28
29 const LO oneLO = Teuchos::OrdinalTraits<LO>::one();
30
31 // Count number of successfully replaced diagonal entries
32 LO numReplacedDiagEntries = 0;
33
34 // Extract some useful data
35 auto rowMapPtr = matrix.getRowMap();
36 if (rowMapPtr.is_null() || rowMapPtr->getComm().is_null()) {
37 // Processes on which the row Map or its communicator is null
38 // don't participate. Users shouldn't even call this method on
39 // those processes.
40 return numReplacedDiagEntries;
41 }
42 auto colMapPtr = matrix.getColMap();
43
44 TEUCHOS_TEST_FOR_EXCEPTION(colMapPtr.get() == nullptr, std::invalid_argument,
45 "replaceDiagonalCrsMatrix: "
46 "Input matrix must have a nonnull column Map.");
47
48 const map_type& rowMap = *rowMapPtr;
49 const map_type& colMap = *colMapPtr;
50 const LO myNumRows = static_cast<LO>(matrix.getLocalNumRows());
51 const bool isFillCompleteOnInput = matrix.isFillComplete();
52
54 TEUCHOS_TEST_FOR_EXCEPTION(!rowMap.isSameAs(*newDiag.getMap()), std::runtime_error,
55 "Row map of matrix and map of input vector do not match.");
56 }
57
58 // KJ: This fence is necessary for UVM. Views used in the row map and colmap
59 // can use UVM and they are accessed in the following routine. So, we need to
60 // make sure that the values are available for touching in host.
61 typename crs_matrix_type::execution_space().fence();
62
63 if (isFillCompleteOnInput)
64 matrix.resumeFill();
65
66 Teuchos::ArrayRCP<const SC> newDiagData = newDiag.getVector(0)->getData();
67 LO numReplacedEntriesPerRow = 0;
68
69 auto invalid = Teuchos::OrdinalTraits<LO>::invalid();
70
71 // Loop over all local rows to replace the diagonal entry row by row
72 for (LO lclRowInd = 0; lclRowInd < myNumRows; ++lclRowInd) {
73 // Get row and column indices of the diagonal entry
74 const GO gblInd = rowMap.getGlobalElement(lclRowInd);
75 const LO lclColInd = colMap.getLocalElement(gblInd);
76
77 // If the row map is not one-to-one, the diagonal may not be on this proc.
78 // Skip this row; some processor will have the diagonal for this row.
79 if (lclColInd == invalid) continue;
80
81 const SC vals[] = {static_cast<SC>(newDiagData[lclRowInd])};
82 const LO cols[] = {lclColInd};
83
84 // Do the actual replacement of the diagonal element, if on this proc
85 numReplacedEntriesPerRow = matrix.replaceLocalValues(lclRowInd, oneLO,
86 vals, cols);
87
88 // Check for success of replacement.
89 // numReplacedEntriesPerRow is one if the diagonal was replaced.
90 // numReplacedEntriesPerRow is zero if the diagonal is not on
91 // this processor. For example, in a 2D matrix distribution, gblInd may
92 // be in both the row and column map, but the diagonal may not be on
93 // this processor.
94 if (numReplacedEntriesPerRow == oneLO) {
95 ++numReplacedDiagEntries;
96 }
97 }
98
99 if (isFillCompleteOnInput)
100 matrix.fillComplete();
101
102 return numReplacedDiagEntries;
103}
104
105} // namespace Tpetra
106//
107// Explicit instantiation macro
108//
109// Must be expanded from within the Tpetra namespace!
110//
111
112#define TPETRA_REPLACEDIAGONALCRSMATRIX_INSTANT(SCALAR, LO, GO, NODE) \
113 \
114 template LO replaceDiagonalCrsMatrix( \
115 CrsMatrix<SCALAR, LO, GO, NODE>& matrix, \
116 const Vector<SCALAR, LO, GO, NODE>& newDiag);
117
118#endif // #ifndef TPETRA_REPLACEDIAGONALCRSMATRIX_DEF_HPP
Declaration of Tpetra::Details::Behavior, a class that describes Tpetra's behavior.
static bool debug()
Whether Tpetra is in debug mode.
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.