Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_DenseContainer_def.hpp
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2009) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 */
42 
43 #ifndef IFPACK2_DENSECONTAINER_DEF_HPP
44 #define IFPACK2_DENSECONTAINER_DEF_HPP
45 
46 #include "Tpetra_CrsMatrix.hpp"
47 #include "Teuchos_LAPACK.hpp"
48 #include "Tpetra_BlockMultiVector.hpp"
49 
50 #ifdef HAVE_MPI
51 # include <mpi.h>
53 #else
54 # include "Teuchos_DefaultSerialComm.hpp"
55 #endif // HAVE_MPI
56 
57 namespace Ifpack2 {
58 
59 template<class MatrixType, class LocalScalarType>
62  const Teuchos::Array<Teuchos::Array<LO> >& partitions,
64  bool pointIndexed) :
65  ContainerImpl<MatrixType, LocalScalarType> (matrix, partitions, pointIndexed),
66  scalarOffsets_ (this->numBlocks_)
67 {
69  !matrix->hasColMap(), std::invalid_argument, "Ifpack2::DenseContainer: "
70  "The constructor's input matrix must have a column Map.");
71 
72  //compute scalarOffsets_
73  GO totalScalars = 0;
74  for(LO i = 0; i < this->numBlocks_; i++)
75  {
76  scalarOffsets_[i] = totalScalars;
77  totalScalars += this->blockSizes_[i] * this->scalarsPerRow_ *
78  this->blockSizes_[i] * this->scalarsPerRow_;
79  }
80  scalars_.resize(totalScalars);
81  for(int i = 0; i < this->numBlocks_; i++)
82  {
83  LO denseRows = this->blockSizes_[i] * this->scalarsPerRow_;
84  //create square dense matrix (stride is same as rows and cols)
85  diagBlocks_.emplace_back(Teuchos::View, scalars_.data() + scalarOffsets_[i], denseRows, denseRows, denseRows);
86  }
87 
88  ipiv_.resize(this->blockRows_.size() * this->scalarsPerRow_);
89 }
90 
91 template<class MatrixType, class LocalScalarType>
94 
95 template<class MatrixType, class LocalScalarType>
96 void
99 {
100  // Fill the diagonal block and LU permutation array with zeros.
101  for(int i = 0; i < this->numBlocks_; i++)
102  diagBlocks_[i].putScalar(Teuchos::ScalarTraits<LSC>::zero());
103  std::fill (ipiv_.begin (), ipiv_.end (), 0);
104 
105  this->IsInitialized_ = true;
106  // We assume that if you called this method, you intend to recompute
107  // everything.
108  this->IsComputed_ = false;
109 }
110 
111 template<class MatrixType, class LocalScalarType>
112 void
115 {
116 // FIXME: I am commenting this out because it breaks block CRS support
117 // TEUCHOS_TEST_FOR_EXCEPTION(
118 // static_cast<size_t> (ipiv_.size ()) != numRows_, std::logic_error,
119 // "Ifpack2::DenseContainer::compute: ipiv_ array has the wrong size. "
120 // "Please report this bug to the Ifpack2 developers.");
121 
122  this->IsComputed_ = false;
123  if (!this->isInitialized ()) {
124  this->initialize();
125  }
126 
127  extract (); // Extract the submatrices
128  factor (); // Factor them
129  this->IsComputed_ = true;
130 }
131 
132 template<class MatrixType, class LocalScalarType>
134 {
135  using Teuchos::Array;
136  using Teuchos::ArrayView;
137  using STS = Teuchos::ScalarTraits<SC>;
138  SC zero = STS::zero();
139  const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
140  //To extract diagonal blocks, need to translate local rows to local columns.
141  //Strategy: make a lookup table that translates local cols in the matrix to offsets in blockRows_:
142  //blockOffsets_[b] <= offset < blockOffsets_[b+1]: tests whether the column is in block b.
143  //offset - blockOffsets_[b]: gives the column within block b.
144  //
145  //This provides the block and col within a block in O(1).
146  if(this->scalarsPerRow_ > 1)
147  {
148  Array<LO> colToBlockOffset(this->inputBlockMatrix_->getNodeNumCols(), INVALID);
149  for(int i = 0; i < this->numBlocks_; i++)
150  {
151  //Get the interval where block i is defined in blockRows_
152  LO blockStart = this->blockOffsets_[i];
153  LO blockEnd = blockStart + this->blockSizes_[i];
154  ArrayView<const LO> blockRows = this->getBlockRows(i);
155  //Set the lookup table entries for the columns appearing in block i.
156  //If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
157  //this is OK. The values updated here are only needed to process block i's entries.
158  for(size_t j = 0; j < size_t(blockRows.size()); j++)
159  {
160  LO localCol = this->translateRowToCol(blockRows[j]);
161  colToBlockOffset[localCol] = blockStart + j;
162  }
163  for(LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++)
164  {
165  //get a raw view of the whole block row
166  const LO* indices;
167  SC* values;
168  LO numEntries;
169  LO inputRow = this->blockRows_[blockStart + blockRow];
170  this->inputBlockMatrix_->getLocalRowView(inputRow, indices, values, numEntries);
171  for(LO k = 0; k < numEntries; k++)
172  {
173  LO colOffset = colToBlockOffset[indices[k]];
174  if(blockStart <= colOffset && colOffset < blockEnd)
175  {
176  //This entry does appear in the diagonal block.
177  //(br, bc) identifies the scalar's position in the BlockCrs block.
178  //Convert this to (r, c) which is its position in the container block.
179  LO blockCol = colOffset - blockStart;
180  for(LO bc = 0; bc < this->bcrsBlockSize_; bc++)
181  {
182  for(LO br = 0; br < this->bcrsBlockSize_; br++)
183  {
184  LO r = this->bcrsBlockSize_ * blockRow + br;
185  LO c = this->bcrsBlockSize_ * blockCol + bc;
186  auto val = values[k * (this->bcrsBlockSize_ * this->bcrsBlockSize_) + (br + this->bcrsBlockSize_ * bc)];
187  if(val != zero)
188  diagBlocks_[i](r, c) = val;
189  }
190  }
191  }
192  }
193  }
194  }
195  }
196  else
197  {
198  //get the mapping from point-indexed matrix columns to offsets in blockRows_
199  //(this includes regular CrsMatrix columns, in which case bcrsBlockSize_ == 1)
200  Array<LO> colToBlockOffset(this->inputMatrix_->getNodeNumCols() * this->bcrsBlockSize_, INVALID);
201  for(int i = 0; i < this->numBlocks_; i++)
202  {
203  //Get the interval where block i is defined in blockRows_
204  LO blockStart = this->blockOffsets_[i];
205  LO blockEnd = blockStart + this->blockSizes_[i];
206  ArrayView<const LO> blockRows = this->getBlockRows(i);
207  //Set the lookup table entries for the columns appearing in block i.
208  //If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
209  //this is OK. The values updated here are only needed to process block i's entries.
210  for(size_t j = 0; j < size_t(blockRows.size()); j++)
211  {
212  //translateRowToCol will return the corresponding split column
213  LO localCol = this->translateRowToCol(blockRows[j]);
214  colToBlockOffset[localCol] = blockStart + j;
215  }
216  for(size_t blockRow = 0; blockRow < size_t(blockRows.size()); blockRow++)
217  {
218  //get a view of the split row
219  LO inputPointRow = this->blockRows_[blockStart + blockRow];
220  auto rowView = this->getInputRowView(inputPointRow);
221  for(size_t k = 0; k < rowView.size(); k++)
222  {
223  LO colOffset = colToBlockOffset[rowView.ind(k)];
224  if(blockStart <= colOffset && colOffset < blockEnd)
225  {
226  LO blockCol = colOffset - blockStart;
227  auto val = rowView.val(k);
228  if(val != zero)
229  diagBlocks_[i](blockRow, blockCol) = rowView.val(k);
230  }
231  }
232  }
233  }
234  }
235 }
236 
237 
238 template<class MatrixType, class LocalScalarType>
239 void
240 DenseContainer<MatrixType, LocalScalarType>::
241 factor ()
242 {
244  for(int i = 0; i < this->numBlocks_; i++)
245  {
246  int INFO = 0;
247  int* blockIpiv = &ipiv_[this->blockOffsets_[i] * this->scalarsPerRow_];
248  lapack.GETRF(diagBlocks_[i].numRows(),
249  diagBlocks_[i].numCols(),
250  diagBlocks_[i].values(),
251  diagBlocks_[i].stride(),
252  blockIpiv, &INFO);
253  // INFO < 0 is a bug.
255  INFO < 0, std::logic_error, "Ifpack2::DenseContainer::factor: "
256  "LAPACK's _GETRF (LU factorization with partial pivoting) was called "
257  "incorrectly. INFO = " << INFO << " < 0. "
258  "Please report this bug to the Ifpack2 developers.");
259  // INFO > 0 means the matrix is singular. This is probably an issue
260  // either with the choice of rows the rows we extracted, or with the
261  // input matrix itself.
263  INFO > 0, std::runtime_error, "Ifpack2::DenseContainer::factor: "
264  "LAPACK's _GETRF (LU factorization with partial pivoting) reports that the "
265  "computed U factor is exactly singular. U(" << INFO << "," << INFO << ") "
266  "(one-based index i) is exactly zero. This probably means that the input "
267  "matrix has a singular diagonal block.");
268  }
269 }
270 
271 template<class MatrixType, class LocalScalarType>
272 void
273 DenseContainer<MatrixType, LocalScalarType>::
274 solveBlock(HostSubviewLocal X,
275  HostSubviewLocal Y,
276  int blockIndex,
277  Teuchos::ETransp mode,
278  LSC alpha,
279  LSC beta) const
280 {
281  #ifdef HAVE_IFPACK2_DEBUG
283  X.extent (0) != Y.extent (0),
284  std::logic_error, "Ifpack2::DenseContainer::solveBlock: X and Y have "
285  "incompatible dimensions (" << X.extent (0) << " resp. "
286  << Y.extent (0) << "). Please report this bug to "
287  "the Ifpack2 developers.");
288 
290  X.extent (1) != Y.extent(1),
291  std::logic_error, "Ifpack2::DenseContainer::solveBlock: X and Y have "
292  "incompatible numbers of vectors (" << X.extent (1) << " resp. "
293  << Y.extent (1) << "). Please report this bug to "
294  "the Ifpack2 developers.");
295  #endif
296 
297  typedef Teuchos::ScalarTraits<LSC> STS;
298  size_t numRows = X.extent(0);
299  size_t numVecs = X.extent(1);
300  if(alpha == STS::zero()) { // don't need to solve the linear system
301  if(beta == STS::zero()) {
302  // Use BLAS AXPY semantics for beta == 0: overwrite, clobbering
303  // any Inf or NaN values in Y (rather than multiplying them by
304  // zero, resulting in NaN values).
305  for(size_t j = 0; j < numVecs; j++)
306  {
307  for(size_t i = 0; i < numRows; i++)
308  Y(i, j) = STS::zero();
309  }
310  }
311  else // beta != 0
312  for(size_t j = 0; j < numVecs; j++)
313  {
314  for(size_t i = 0; i < numRows; i++)
315  Y(i, j) = beta * Y(i, j);
316  }
317  }
318  else { // alpha != 0; must solve the linear system
320  // If beta is nonzero or Y is not constant stride, we have to use
321  // a temporary output multivector. It gets a (deep) copy of X,
322  // since GETRS overwrites its (multi)vector input with its output.
323  std::vector<LSC> yTemp(numVecs * numRows);
324  for(size_t j = 0; j < numVecs; j++)
325  {
326  for(size_t i = 0; i < numRows; i++)
327  yTemp[j * numRows + i] = X(i, j);
328  }
329 
330  int INFO = 0;
331  int* blockIpiv = &ipiv_[this->blockOffsets_[blockIndex] * this->scalarsPerRow_];
332  const char trans =
333  (mode == Teuchos::CONJ_TRANS ? 'C' : (mode == Teuchos::TRANS ? 'T' : 'N'));
334  lapack.GETRS (trans,
335  diagBlocks_[blockIndex].numRows (),
336  numVecs,
337  diagBlocks_[blockIndex].values (),
338  diagBlocks_[blockIndex].stride (),
339  blockIpiv,
340  yTemp.data(),
341  numRows,
342  &INFO);
343 
344  if (beta != STS::zero ()) {
345  for(size_t j = 0; j < numVecs; j++)
346  {
347  for(size_t i = 0; i < numRows; i++)
348  {
349  Y(i, j) *= ISC(beta);
350  Y(i, j) += ISC(alpha * yTemp[j * numRows + i]);
351  }
352  }
353  }
354  else {
355  for(size_t j = 0; j < numVecs; j++)
356  {
357  for(size_t i = 0; i < numRows; i++)
358  Y(i, j) = ISC(yTemp[j * numRows + i]);
359  }
360  }
361 
363  INFO != 0, std::runtime_error, "Ifpack2::DenseContainer::solveBlock: "
364  "LAPACK's _GETRS (solve using LU factorization with partial pivoting) "
365  "failed with INFO = " << INFO << " != 0.");
366  }
367 }
368 
369 template<class MatrixType, class LocalScalarType>
370 std::ostream&
372 print (std::ostream& os) const
373 {
374  Teuchos::FancyOStream fos (Teuchos::rcpFromRef (os));
375  fos.setOutputToRootOnly (0);
376  this->describe (fos);
377  return os;
378 }
379 
380 template<class MatrixType, class LocalScalarType>
381 std::string
383 description () const
384 {
385  std::ostringstream oss;
386  oss << "Ifpack::DenseContainer: ";
387  if (this->isInitialized()) {
388  if (this->isComputed()) {
389  oss << "{status = initialized, computed";
390  }
391  else {
392  oss << "{status = initialized, not computed";
393  }
394  }
395  else {
396  oss << "{status = not initialized, not computed";
397  }
398 
399  oss << "}";
400  return oss.str();
401 }
402 
403 template<class MatrixType, class LocalScalarType>
404 void
407  const Teuchos::EVerbosityLevel verbLevel) const
408 {
409  using std::endl;
410  if(verbLevel==Teuchos::VERB_NONE) return;
411  os << "================================================================================" << endl;
412  os << "Ifpack2::DenseContainer" << endl;
413  for(int i = 0; i < this->numBlocks_; i++)
414  {
415  os << "Block " << i << " number of rows = " << this->blockSizes_[i] << endl;
416  }
417  os << "isInitialized() = " << this->IsInitialized_ << endl;
418  os << "isComputed() = " << this->IsComputed_ << endl;
419  os << "================================================================================" << endl;
420  os << endl;
421 }
422 
423 template<class MatrixType, class LocalScalarType>
425 {
426  diagBlocks_.clear();
427  scalars_.clear();
429 }
430 
431 template<class MatrixType, class LocalScalarType>
433 {
434  return "Dense";
435 }
436 
437 } // namespace Ifpack2
438 
439 // There's no need to instantiate for CrsMatrix too. All Ifpack2
440 // preconditioners can and should do dynamic casts if they need a type
441 // more specific than RowMatrix.
442 
443 #define IFPACK2_DENSECONTAINER_INSTANT(S,LO,GO,N) \
444  template class Ifpack2::DenseContainer< Tpetra::RowMatrix<S, LO, GO, N>, S >;
445 
446 #endif // IFPACK2_DENSECONTAINER_DEF_HPP
Teuchos::Array< LO > blockRows_
Local indices of the rows of the input matrix that belong to this block.
Definition: Ifpack2_Container_decl.hpp:293
Teuchos::Array< LO > blockSizes_
Number of rows in each block.
Definition: Ifpack2_Container_decl.hpp:295
int numBlocks_
The number of blocks (partitions) in the container.
Definition: Ifpack2_Container_decl.hpp:291
The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS). This class allows a custom scalar type (LocalScalarType) to be used for storing blocks and solving the block systems. Hiding this template parameter from the Container interface simplifies the BlockRelaxation and ContainerFactory classes.
Definition: Ifpack2_Container_decl.hpp:342
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void compute()
Extract the local diagonal block and prepare the solver.
Definition: Ifpack2_DenseContainer_def.hpp:114
LO scalarsPerRow_
Definition: Ifpack2_Container_decl.hpp:316
static std::string getName()
Get the name of this container type for Details::constructContainer()
Definition: Ifpack2_DenseContainer_def.hpp:432
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:61
void GETRF(const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, OrdinalType *info) const
Store and solve a local dense linear problem.
Definition: Ifpack2_DenseContainer_decl.hpp:103
virtual std::ostream & print(std::ostream &os) const
Print information about this object to the given output stream.
Definition: Ifpack2_DenseContainer_def.hpp:372
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:406
basic_FancyOStream & setOutputToRootOnly(const int rootRank)
void resize(size_type new_size, const value_type &x=value_type())
virtual ~DenseContainer()
Destructor (declared virtual for memory safety of derived classes).
Definition: Ifpack2_DenseContainer_def.hpp:93
void GETRS(const char &TRANS, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
size_type size() const
virtual std::string description() const
A one-line description of this object.
Definition: Ifpack2_DenseContainer_def.hpp:383
virtual void initialize()
Do all set-up operations that only require matrix structure.
Definition: Ifpack2_DenseContainer_def.hpp:98