Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_DenseContainer_def.hpp
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef IFPACK2_DENSECONTAINER_DEF_HPP
11#define IFPACK2_DENSECONTAINER_DEF_HPP
12
13#include "Tpetra_CrsMatrix.hpp"
14#include "Teuchos_LAPACK.hpp"
15#include "Tpetra_BlockMultiVector.hpp"
17
18#ifdef HAVE_MPI
19#include <mpi.h>
20#include "Teuchos_DefaultMpiComm.hpp"
21#else
22#include "Teuchos_DefaultSerialComm.hpp"
23#endif // HAVE_MPI
24
25namespace Ifpack2 {
26
27template <class MatrixType, class LocalScalarType>
29 DenseContainer(const Teuchos::RCP<const row_matrix_type>& matrix,
30 const Teuchos::Array<Teuchos::Array<LO> >& partitions,
31 const Teuchos::RCP<const import_type>&,
32 bool pointIndexed)
33 : ContainerImpl<MatrixType, LocalScalarType>(matrix, partitions, pointIndexed)
34 , scalarOffsets_(this->numBlocks_) {
35 TEUCHOS_TEST_FOR_EXCEPTION(
36 !matrix->hasColMap(), std::invalid_argument,
37 "Ifpack2::DenseContainer: "
38 "The constructor's input matrix must have a column Map.");
39
40 // compute scalarOffsets_
41 GO totalScalars = 0;
42 for (LO i = 0; i < this->numBlocks_; i++) {
43 scalarOffsets_[i] = totalScalars;
44 totalScalars += this->blockSizes_[i] * this->scalarsPerRow_ *
45 this->blockSizes_[i] * this->scalarsPerRow_;
46 }
47 scalars_.resize(totalScalars);
48 for (int i = 0; i < this->numBlocks_; i++) {
49 LO denseRows = this->blockSizes_[i] * this->scalarsPerRow_;
50 // create square dense matrix (stride is same as rows and cols)
51 diagBlocks_.emplace_back(Teuchos::View, scalars_.data() + scalarOffsets_[i], denseRows, denseRows, denseRows);
52 }
53
54 ipiv_.resize(this->blockRows_.size() * this->scalarsPerRow_);
55}
56
57template <class MatrixType, class LocalScalarType>
60
61template <class MatrixType, class LocalScalarType>
64 // Fill the diagonal block and LU permutation array with zeros.
65 for (int i = 0; i < this->numBlocks_; i++)
66 diagBlocks_[i].putScalar(Teuchos::ScalarTraits<LSC>::zero());
67 std::fill(ipiv_.begin(), ipiv_.end(), 0);
68
69 this->IsInitialized_ = true;
70 // We assume that if you called this method, you intend to recompute
71 // everything.
72 this->IsComputed_ = false;
73}
74
75template <class MatrixType, class LocalScalarType>
77 compute() {
78 // FIXME: I am commenting this out because it breaks block CRS support
79 // TEUCHOS_TEST_FOR_EXCEPTION(
80 // static_cast<size_t> (ipiv_.size ()) != numRows_, std::logic_error,
81 // "Ifpack2::DenseContainer::compute: ipiv_ array has the wrong size. "
82 // "Please report this bug to the Ifpack2 developers.");
83
84 this->IsComputed_ = false;
85 if (!this->isInitialized()) {
86 this->initialize();
87 }
88
89 extract(); // Extract the submatrices
90 factor(); // Factor them
91 this->IsComputed_ = true;
92}
93
94template <class MatrixType, class LocalScalarType>
96 using Teuchos::Array;
97 using Teuchos::ArrayView;
98 using STS = Teuchos::ScalarTraits<SC>;
99 SC zero = STS::zero();
100 const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
101 // To extract diagonal blocks, need to translate local rows to local columns.
102 // Strategy: make a lookup table that translates local cols in the matrix to offsets in blockRows_:
103 // blockOffsets_[b] <= offset < blockOffsets_[b+1]: tests whether the column is in block b.
104 // offset - blockOffsets_[b]: gives the column within block b.
105 //
106 // This provides the block and col within a block in O(1).
107 if (this->scalarsPerRow_ > 1) {
108 Array<LO> colToBlockOffset(this->inputBlockMatrix_->getLocalNumCols(), INVALID);
109 for (int i = 0; i < this->numBlocks_; i++) {
110 // Get the interval where block i is defined in blockRows_
111 LO blockStart = this->blockOffsets_[i];
112 LO blockEnd = blockStart + this->blockSizes_[i];
113 ArrayView<const LO> blockRows = this->getBlockRows(i);
114 // Set the lookup table entries for the columns appearing in block i.
115 // If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
116 // this is OK. The values updated here are only needed to process block i's entries.
117 for (size_t j = 0; j < size_t(blockRows.size()); j++) {
118 LO localCol = this->translateRowToCol(blockRows[j]);
119 colToBlockOffset[localCol] = blockStart + j;
120 }
121 using h_inds_type = typename block_crs_matrix_type::local_inds_host_view_type;
122 using h_vals_type = typename block_crs_matrix_type::values_host_view_type;
123 for (LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++) {
124 // get a raw view of the whole block row
125 h_inds_type indices;
126 h_vals_type values;
127 LO inputRow = this->blockRows_[blockStart + blockRow];
128 this->inputBlockMatrix_->getLocalRowView(inputRow, indices, values);
129 LO numEntries = (LO)indices.size();
130 for (LO k = 0; k < numEntries; k++) {
131 LO colOffset = colToBlockOffset[indices[k]];
132 if (blockStart <= colOffset && colOffset < blockEnd) {
133 // This entry does appear in the diagonal block.
134 //(br, bc) identifies the scalar's position in the BlockCrs block.
135 // Convert this to (r, c) which is its position in the container block.
136 LO blockCol = colOffset - blockStart;
137 for (LO bc = 0; bc < this->bcrsBlockSize_; bc++) {
138 for (LO br = 0; br < this->bcrsBlockSize_; br++) {
139 LO r = this->bcrsBlockSize_ * blockRow + br;
140 LO c = this->bcrsBlockSize_ * blockCol + bc;
141 auto val = values[k * (this->bcrsBlockSize_ * this->bcrsBlockSize_) + (br + this->bcrsBlockSize_ * bc)];
142 if (val != zero)
143 diagBlocks_[i](r, c) = val;
144 }
145 }
146 }
147 }
148 }
149 }
150 } else {
151 // get the mapping from point-indexed matrix columns to offsets in blockRows_
152 //(this includes regular CrsMatrix columns, in which case bcrsBlockSize_ == 1)
153 Array<LO> colToBlockOffset(this->inputMatrix_->getLocalNumCols() * this->bcrsBlockSize_, INVALID);
154 for (int i = 0; i < this->numBlocks_; i++) {
155 // Get the interval where block i is defined in blockRows_
156 LO blockStart = this->blockOffsets_[i];
157 LO blockEnd = blockStart + this->blockSizes_[i];
158 ArrayView<const LO> blockRows = this->getBlockRows(i);
159 // Set the lookup table entries for the columns appearing in block i.
160 // If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
161 // this is OK. The values updated here are only needed to process block i's entries.
162 for (size_t j = 0; j < size_t(blockRows.size()); j++) {
163 // translateRowToCol will return the corresponding split column
164 LO localCol = this->translateRowToCol(blockRows[j]);
165 colToBlockOffset[localCol] = blockStart + j;
166 }
167 for (size_t blockRow = 0; blockRow < size_t(blockRows.size()); blockRow++) {
168 // get a view of the split row
169 LO inputPointRow = this->blockRows_[blockStart + blockRow];
170 auto rowView = this->getInputRowView(inputPointRow);
171 for (size_t k = 0; k < rowView.size(); k++) {
172 LO colOffset = colToBlockOffset[rowView.ind(k)];
173 if (blockStart <= colOffset && colOffset < blockEnd) {
174 LO blockCol = colOffset - blockStart;
175 auto val = rowView.val(k);
176 if (val != zero)
177 diagBlocks_[i](blockRow, blockCol) = rowView.val(k);
178 }
179 }
180 }
181 }
182 }
183}
184
185template <class MatrixType, class LocalScalarType>
186void DenseContainer<MatrixType, LocalScalarType>::
187 factor() {
188 Teuchos::LAPACK<int, LSC> lapack;
189 for (int i = 0; i < this->numBlocks_; i++) {
190 int INFO = 0;
191 int* blockIpiv = &ipiv_[this->blockOffsets_[i] * this->scalarsPerRow_];
192 lapack.GETRF(diagBlocks_[i].numRows(),
193 diagBlocks_[i].numCols(),
194 diagBlocks_[i].values(),
195 diagBlocks_[i].stride(),
196 blockIpiv, &INFO);
197 // INFO < 0 is a bug.
198 TEUCHOS_TEST_FOR_EXCEPTION(
199 INFO < 0, std::logic_error,
200 "Ifpack2::DenseContainer::factor: "
201 "LAPACK's _GETRF (LU factorization with partial pivoting) was called "
202 "incorrectly. INFO = "
203 << INFO << " < 0. "
204 "Please report this bug to the Ifpack2 developers.");
205 // INFO > 0 means the matrix is singular. This is probably an issue
206 // either with the choice of rows the rows we extracted, or with the
207 // input matrix itself.
208 TEUCHOS_TEST_FOR_EXCEPTION(
209 INFO > 0, std::runtime_error,
210 "Ifpack2::DenseContainer::factor: "
211 "LAPACK's _GETRF (LU factorization with partial pivoting) reports that the "
212 "computed U factor is exactly singular. U("
213 << INFO << "," << INFO << ") "
214 "(one-based index i) is exactly zero. This probably means that the input "
215 "matrix has a singular diagonal block.");
216 }
217}
218
219template <class MatrixType, class LocalScalarType>
220void DenseContainer<MatrixType, LocalScalarType>::
221 solveBlock(ConstHostSubviewLocal X,
222 HostSubviewLocal Y,
223 int blockIndex,
224 Teuchos::ETransp mode,
225 LSC alpha,
226 LSC beta) const {
228 TEUCHOS_TEST_FOR_EXCEPTION(
229 X.extent(0) != Y.extent(0),
230 std::logic_error,
231 "Ifpack2::DenseContainer::solveBlock: X and Y have "
232 "incompatible dimensions ("
233 << X.extent(0) << " resp. "
234 << Y.extent(0) << "). Please report this bug to "
235 "the Ifpack2 developers.");
236
237 TEUCHOS_TEST_FOR_EXCEPTION(
238 X.extent(1) != Y.extent(1),
239 std::logic_error,
240 "Ifpack2::DenseContainer::solveBlock: X and Y have "
241 "incompatible numbers of vectors ("
242 << X.extent(1) << " resp. "
243 << Y.extent(1) << "). Please report this bug to "
244 "the Ifpack2 developers.");
245 }
246
247 typedef Teuchos::ScalarTraits<LSC> STS;
248 size_t numRows = X.extent(0);
249 size_t numVecs = X.extent(1);
250 if (alpha == STS::zero()) { // don't need to solve the linear system
251 if (beta == STS::zero()) {
252 // Use BLAS AXPY semantics for beta == 0: overwrite, clobbering
253 // any Inf or NaN values in Y (rather than multiplying them by
254 // zero, resulting in NaN values).
255 for (size_t j = 0; j < numVecs; j++) {
256 for (size_t i = 0; i < numRows; i++)
257 Y(i, j) = STS::zero();
258 }
259 } else // beta != 0
260 for (size_t j = 0; j < numVecs; j++) {
261 for (size_t i = 0; i < numRows; i++)
262 Y(i, j) = beta * Y(i, j);
263 }
264 } else { // alpha != 0; must solve the linear system
265 Teuchos::LAPACK<int, LSC> lapack;
266 // If beta is nonzero or Y is not constant stride, we have to use
267 // a temporary output multivector. It gets a (deep) copy of X,
268 // since GETRS overwrites its (multi)vector input with its output.
269 std::vector<LSC> yTemp(numVecs * numRows);
270 for (size_t j = 0; j < numVecs; j++) {
271 for (size_t i = 0; i < numRows; i++)
272 yTemp[j * numRows + i] = X(i, j);
273 }
274
275 int INFO = 0;
276 int* blockIpiv = &ipiv_[this->blockOffsets_[blockIndex] * this->scalarsPerRow_];
277 const char trans =
278 (mode == Teuchos::CONJ_TRANS ? 'C' : (mode == Teuchos::TRANS ? 'T' : 'N'));
279 lapack.GETRS(trans,
280 diagBlocks_[blockIndex].numRows(),
281 numVecs,
282 diagBlocks_[blockIndex].values(),
283 diagBlocks_[blockIndex].stride(),
284 blockIpiv,
285 yTemp.data(),
286 numRows,
287 &INFO);
288
289 if (beta != STS::zero()) {
290 for (size_t j = 0; j < numVecs; j++) {
291 for (size_t i = 0; i < numRows; i++) {
292 Y(i, j) *= ISC(beta);
293 Y(i, j) += ISC(alpha * yTemp[j * numRows + i]);
294 }
295 }
296 } else {
297 for (size_t j = 0; j < numVecs; j++) {
298 for (size_t i = 0; i < numRows; i++)
299 Y(i, j) = ISC(yTemp[j * numRows + i]);
300 }
301 }
302
303 TEUCHOS_TEST_FOR_EXCEPTION(
304 INFO != 0, std::runtime_error,
305 "Ifpack2::DenseContainer::solveBlock: "
306 "LAPACK's _GETRS (solve using LU factorization with partial pivoting) "
307 "failed with INFO = "
308 << INFO << " != 0.");
309 }
310}
311
312template <class MatrixType, class LocalScalarType>
313std::ostream&
315 print(std::ostream& os) const {
316 Teuchos::FancyOStream fos(Teuchos::rcpFromRef(os));
317 fos.setOutputToRootOnly(0);
318 this->describe(fos);
319 return os;
320}
321
322template <class MatrixType, class LocalScalarType>
323std::string
325 description() const {
326 std::ostringstream oss;
327 oss << "Ifpack::DenseContainer: ";
328 if (this->isInitialized()) {
329 if (this->isComputed()) {
330 oss << "{status = initialized, computed";
331 } else {
332 oss << "{status = initialized, not computed";
333 }
334 } else {
335 oss << "{status = not initialized, not computed";
336 }
337
338 oss << "}";
339 return oss.str();
340}
341
342template <class MatrixType, class LocalScalarType>
344 describe(Teuchos::FancyOStream& os,
345 const Teuchos::EVerbosityLevel verbLevel) const {
346 using std::endl;
347 if (verbLevel == Teuchos::VERB_NONE) return;
348 os << "================================================================================" << endl;
349 os << "Ifpack2::DenseContainer" << endl;
350 for (int i = 0; i < this->numBlocks_; i++) {
351 os << "Block " << i << " number of rows = " << this->blockSizes_[i] << endl;
352 }
353 os << "isInitialized() = " << this->IsInitialized_ << endl;
354 os << "isComputed() = " << this->IsComputed_ << endl;
355 os << "================================================================================" << endl;
356 os << endl;
357}
358
359template <class MatrixType, class LocalScalarType>
361 diagBlocks_.clear();
362 scalars_.clear();
364}
365
366template <class MatrixType, class LocalScalarType>
368 return "Dense";
369}
370
371} // namespace Ifpack2
372
373// There's no need to instantiate for CrsMatrix too. All Ifpack2
374// preconditioners can and should do dynamic casts if they need a type
375// more specific than RowMatrix.
376
377#define IFPACK2_DENSECONTAINER_INSTANT(S, LO, GO, N) \
378 template class Ifpack2::DenseContainer<Tpetra::RowMatrix<S, LO, GO, N>, S>;
379
380#endif // IFPACK2_DENSECONTAINER_DEF_HPP
Declaration of Ifpack2::Details::Behavior, a class that describes Ifpack2's run-time behavior.
int numBlocks_
The number of blocks (partitions) in the container.
Definition Ifpack2_Container_decl.hpp:257
LO scalarsPerRow_
Definition Ifpack2_Container_decl.hpp:282
Teuchos::Array< LO > blockSizes_
Number of rows in each block.
Definition Ifpack2_Container_decl.hpp:261
Teuchos::Array< LO > blockRows_
Local indices of the rows of the input matrix that belong to this block.
Definition Ifpack2_Container_decl.hpp:259
The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS)....
Definition Ifpack2_Container_decl.hpp:307
Store and solve a local dense linear problem.
Definition Ifpack2_DenseContainer_decl.hpp:71
virtual std::string description() const
A one-line description of this object.
Definition Ifpack2_DenseContainer_def.hpp:325
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to the given FancyOStream.
Definition Ifpack2_DenseContainer_def.hpp:344
static std::string getName()
Get the name of this container type for Details::constructContainer()
Definition Ifpack2_DenseContainer_def.hpp:367
virtual std::ostream & print(std::ostream &os) const
Print information about this object to the given output stream.
Definition Ifpack2_DenseContainer_def.hpp:315
virtual void compute()
Extract the local diagonal block and prepare the solver.
Definition Ifpack2_DenseContainer_def.hpp:77
DenseContainer(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO > > &partitions, const Teuchos::RCP< const import_type > &importer, bool pointIndexed)
Constructor.
Definition Ifpack2_DenseContainer_def.hpp:29
virtual ~DenseContainer()
Destructor (declared virtual for memory safety of derived classes).
Definition Ifpack2_DenseContainer_def.hpp:59
virtual void initialize()
Do all set-up operations that only require matrix structure.
Definition Ifpack2_DenseContainer_def.hpp:63
static bool debug()
Whether Ifpack2 is in debug mode.
Definition Ifpack2_Details_Behavior.cpp:31
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40