10#ifndef IFPACK2_DETAILS_DENSESOLVER_DEF_HPP
11#define IFPACK2_DETAILS_DENSESOLVER_DEF_HPP
13#include "Ifpack2_LocalFilter.hpp"
14#include "Teuchos_LAPACK.hpp"
15#include "Ifpack2_Details_DenseSolver.hpp"
16#include "Tpetra_Map.hpp"
20#include "Teuchos_DefaultMpiComm.hpp"
22#include "Teuchos_DefaultSerialComm.hpp"
32template <
class MatrixType>
34 DenseSolver(
const Teuchos::RCP<const row_matrix_type>& A)
36 , initializeTime_(0.0)
42 , isInitialized_(
false)
43 , isComputed_(
false) {}
45template <
class MatrixType>
46Teuchos::RCP<const typename DenseSolver<MatrixType, false>::map_type>
49 A_.is_null(), std::runtime_error,
50 "Ifpack2::Details::DenseSolver::"
51 "getDomainMap: The input matrix A is null. Please call setMatrix() with a "
52 "nonnull input matrix before calling this method.");
55 return A_->getRangeMap();
58template <
class MatrixType>
59Teuchos::RCP<const typename DenseSolver<MatrixType, false>::map_type>
62 A_.is_null(), std::runtime_error,
63 "Ifpack2::Details::DenseSolver::"
64 "getRangeMap: The input matrix A is null. Please call setMatrix() with a "
65 "nonnull input matrix before calling this method.");
68 return A_->getDomainMap();
71template <
class MatrixType>
77template <
class MatrixType>
79 return isInitialized_;
82template <
class MatrixType>
87template <
class MatrixType>
89 return numInitialize_;
92template <
class MatrixType>
97template <
class MatrixType>
102template <
class MatrixType>
105 return initializeTime_;
108template <
class MatrixType>
114template <
class MatrixType>
120template <
class MatrixType>
121Teuchos::RCP<const typename DenseSolver<MatrixType, false>::row_matrix_type>
126template <
class MatrixType>
129 isInitialized_ =
false;
131 A_local_ = Teuchos::null;
132 A_local_dense_.reshape(0, 0);
136template <
class MatrixType>
138 setMatrix(
const Teuchos::RCP<const row_matrix_type>& A) {
142 const global_size_t numRows = A->getRangeMap()->getGlobalNumElements();
143 const global_size_t
numCols = A->getDomainMap()->getGlobalNumElements();
145 numRows !=
numCols, std::invalid_argument,
146 "Ifpack2::Details::DenseSolver::"
147 "setMatrix: Input matrix must be (globally) square. "
148 "The matrix you provided is "
149 << numRows <<
" by " <<
numCols <<
".");
158template <
class MatrixType>
165 using Teuchos::TimeMonitor;
166 const std::string
timerName(
"Ifpack2::Details::DenseSolver::initialize");
169 if (
timer.is_null()) {
179 A_.is_null(), std::runtime_error,
180 "Ifpack2::Details::DenseSolver::"
181 "initialize: The input matrix A is null. Please call setMatrix() "
182 "with a nonnull input before calling this method.");
185 !A_->hasColMap(), std::invalid_argument,
186 "Ifpack2::Details::DenseSolver: "
187 "The constructor's input matrix must have a column Map, "
188 "so that it has local indices.");
194 if (A_->getComm()->getSize() > 1) {
201 A_local_.is_null(), std::logic_error,
202 "Ifpack2::Details::DenseSolver::"
203 "initialize: A_local_ is null after it was supposed to have been "
204 "initialized. Please report this bug to the Ifpack2 developers.");
207 const size_t numRows = A_local_->getLocalNumRows();
208 const size_t numCols = A_local_->getLocalNumCols();
210 numRows !=
numCols, std::logic_error,
211 "Ifpack2::Details::DenseSolver::"
212 "initialize: Local filter matrix is not square. This should never happen. "
213 "Please report this bug to the Ifpack2 developers.");
214 A_local_dense_.reshape(numRows,
numCols);
215 ipiv_.resize(std::min(numRows,
numCols));
216 std::fill(ipiv_.begin(), ipiv_.end(), 0);
218 isInitialized_ =
true;
225template <
class MatrixType>
228template <
class MatrixType>
231 const std::string
timerName(
"Ifpack2::Details::DenseSolver::compute");
234 if (
timer.is_null()) {
244 A_.is_null(), std::runtime_error,
245 "Ifpack2::Details::DenseSolver::"
246 "compute: The input matrix A is null. Please call setMatrix() with a "
247 "nonnull input, then call initialize(), before calling this method.");
250 A_local_.is_null(), std::logic_error,
251 "Ifpack2::Details::DenseSolver::"
252 "compute: A_local_ is null. Please report this bug to the Ifpack2 "
256 if (!this->isInitialized()) {
259 extract(A_local_dense_, *A_local_);
260 factor(A_local_dense_, ipiv_());
268template <
class MatrixType>
270 factor(Teuchos::SerialDenseMatrix<int, scalar_type>& A,
271 const Teuchos::ArrayView<int>& ipiv) {
273 std::fill(ipiv.begin(), ipiv.end(), 0);
275 Teuchos::LAPACK<int, scalar_type>
lapack;
277 lapack.GETRF(A.numRows(), A.numCols(), A.values(), A.stride(),
278 ipiv.getRawPtr(), &
INFO);
281 INFO < 0, std::logic_error,
282 "Ifpack2::Details::DenseSolver::factor: "
283 "LAPACK's _GETRF (LU factorization with partial pivoting) was called "
284 "incorrectly. INFO = "
286 "Please report this bug to the Ifpack2 developers.");
291 INFO > 0, std::runtime_error,
292 "Ifpack2::Details::DenseSolver::factor: "
293 "LAPACK's _GETRF (LU factorization with partial pivoting) reports that the "
294 "computed U factor is exactly singular. U("
296 "(one-based index i) is exactly zero. This probably means that the input "
297 "matrix has a singular diagonal block.");
300template <
class MatrixType>
301void DenseSolver<MatrixType, false>::
302 applyImpl(
const MV& X,
304 const Teuchos::ETransp mode,
305 const scalar_type alpha,
306 const scalar_type beta)
const {
307 using Teuchos::ArrayRCP;
308 using Teuchos::CONJ_TRANS;
311 using Teuchos::rcpFromRef;
312 using Teuchos::TRANS;
314 const int numVecs =
static_cast<int>(X.getNumVectors());
315 if (alpha == STS::zero()) {
316 if (beta == STS::zero()) {
320 Y.putScalar(STS::zero());
322 Y.scale(STS::zero());
325 Teuchos::LAPACK<int, scalar_type> lapack;
331 if (beta == STS::zero() && Y.isConstantStride() && alpha == STS::one()) {
333 Y_tmp = rcpFromRef(Y);
335 Y_tmp = rcp(
new MV(X, Teuchos::Copy));
336 if (alpha != STS::one()) {
340 const int Y_stride =
static_cast<int>(Y_tmp->getStride());
341 ArrayRCP<scalar_type> Y_view = Y_tmp->get1dViewNonConst();
342 scalar_type*
const Y_ptr = Y_view.getRawPtr();
345 (mode == CONJ_TRANS ?
'C' : (mode ==
TRANS ?
'T' :
'N'));
346 lapack.GETRS(trans, A_local_dense_.numRows(), numVecs,
347 A_local_dense_.values(), A_local_dense_.stride(),
348 ipiv_.getRawPtr(), Y_ptr, Y_stride, &INFO);
349 TEUCHOS_TEST_FOR_EXCEPTION(
350 INFO != 0, std::runtime_error,
351 "Ifpack2::Details::DenseSolver::"
352 "applyImpl: LAPACK's _GETRS (solve using LU factorization with "
353 "partial pivoting) failed with INFO = "
354 << INFO <<
" != 0.");
356 if (beta != STS::zero()) {
357 Y.update(alpha, *Y_tmp, beta);
358 }
else if (!Y.isConstantStride()) {
359 deep_copy(Y, *Y_tmp);
364template <
class MatrixType>
366 apply(
const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
367 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
368 Teuchos::ETransp
mode,
371 using Teuchos::ArrayView;
375 using Teuchos::rcpFromRef;
377 const std::string
timerName(
"Ifpack2::Details::DenseSolver::apply");
379 if (
timer.is_null()) {
390 !isComputed_, std::runtime_error,
391 "Ifpack2::Details::DenseSolver::apply: "
392 "You must have called the compute() method before you may call apply(). "
393 "You may call the apply() method as many times as you want after calling "
394 "compute() once, but you must have called compute() at least once.");
396 const size_t numVecs = X.getNumVectors();
399 numVecs != Y.getNumVectors(), std::runtime_error,
400 "Ifpack2::Details::DenseSolver::apply: X and Y have different numbers "
402 << X.getNumVectors() <<
", but Y has "
403 << X.getNumVectors() <<
".");
412 const bool multipleProcs = (A_->getRowMap()->getComm()->getSize() >= 1);
417 X_local = X.offsetView(A_local_->getDomainMap(), 0);
418 Y_local = Y.offsetViewNonConst(A_local_->getRangeMap(), 0);
435template <
class MatrixType>
438 std::ostringstream
out;
443 out <<
"\"Ifpack2::Details::DenseSolver\": ";
448 out <<
"Initialized: " << (isInitialized() ?
"true" :
"false") <<
", "
449 <<
"Computed: " << (isComputed() ?
"true" :
"false") <<
", ";
452 out <<
"Matrix: null";
454 out <<
"Matrix: not null"
455 <<
", Global matrix dimensions: ["
456 << A_->getGlobalNumRows() <<
", " << A_->getGlobalNumCols() <<
"]";
463template <
class MatrixType>
466 const Teuchos::EVerbosityLevel
verbLevel)
const {
468 using Teuchos::FancyOStream;
469 using Teuchos::OSTab;
471 using Teuchos::rcpFromRef;
476 RCP<FancyOStream> ptrOut = rcpFromRef(out);
478 if (this->getObjectLabel() !=
"") {
479 out <<
"label: " << this->getObjectLabel() << endl;
481 out <<
"initialized: " << (isInitialized_ ?
"true" :
"false") << endl
482 <<
"computed: " << (isComputed_ ?
"true" :
"false") << endl
483 <<
"number of initialize calls: " << numInitialize_ << endl
484 <<
"number of compute calls: " << numCompute_ << endl
485 <<
"number of apply calls: " << numApply_ << endl
486 <<
"total time in seconds in initialize: " << initializeTime_ << endl
487 <<
"total time in seconds in compute: " << computeTime_ << endl
488 <<
"total time in seconds in apply: " << applyTime_ << endl;
489 if (verbLevel >= Teuchos::VERB_EXTREME) {
490 out <<
"A_local_dense_:" << endl;
494 for (
int i = 0; i < A_local_dense_.numRows(); ++i) {
495 for (
int j = 0; j < A_local_dense_.numCols(); ++j) {
496 out << A_local_dense_(i, j);
497 if (j + 1 < A_local_dense_.numCols()) {
501 if (i + 1 < A_local_dense_.numRows()) {
507 out <<
"ipiv_: " << Teuchos::toString(ipiv_) << endl;
512template <
class MatrixType>
515 const Teuchos::EVerbosityLevel
verbLevel)
const {
517 using Teuchos::FancyOStream;
518 using Teuchos::OSTab;
520 using Teuchos::rcpFromRef;
530 out <<
"Ifpack2::Details::DenseSolver:" <<
endl;
536 const Teuchos::Comm<int>&
comm = *(A_->getRowMap()->getComm());
540 out <<
"Ifpack2::Details::DenseSolver:" <<
endl;
555template <
class MatrixType>
558 const row_matrix_type&
A_local) {
559 using Teuchos::Array;
560 using Teuchos::ArrayView;
561 typedef local_ordinal_type LO;
562 typedef typename Teuchos::ArrayView<LO>::size_type size_type;
578 static_cast<size_type
>(
A_local.getLocalMaxNumRowEntries());
603 const scalar_type val = values[
k];
615template <
class MatrixType>
617 DenseSolver(
const Teuchos::RCP<const row_matrix_type>& A) {
621template <
class MatrixType>
622Teuchos::RCP<const typename DenseSolver<MatrixType, true>::map_type>
627template <
class MatrixType>
628Teuchos::RCP<const typename DenseSolver<MatrixType, true>::map_type>
633template <
class MatrixType>
639template <
class MatrixType>
644template <
class MatrixType>
649template <
class MatrixType>
654template <
class MatrixType>
659template <
class MatrixType>
664template <
class MatrixType>
670template <
class MatrixType>
676template <
class MatrixType>
682template <
class MatrixType>
683Teuchos::RCP<const typename DenseSolver<MatrixType, true>::row_matrix_type>
688template <
class MatrixType>
690 setMatrix(
const Teuchos::RCP<const row_matrix_type>& A) {
694template <
class MatrixType>
699template <
class MatrixType>
704template <
class MatrixType>
709template <
class MatrixType>
711 apply(
const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
712 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
713 Teuchos::ETransp
mode,
719template <
class MatrixType>
725template <
class MatrixType>
728 const Teuchos::EVerbosityLevel
verbLevel)
const {
735#define IFPACK2_DETAILS_DENSESOLVER_INSTANT(S, LO, GO, N) \
736 template class Ifpack2::Details::DenseSolver<Tpetra::RowMatrix<S, LO, GO, N> >;
virtual void setMatrix(const Teuchos::RCP< const Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > &A)=0
Set the new matrix.
MatrixType::scalar_type scalar_type
The type of entries in the input (global) matrix.
Definition Ifpack2_Details_DenseSolver_decl.hpp:71
MatrixType::scalar_type scalar_type
The type of entries in the input (global) matrix.
Definition Ifpack2_Details_DenseSolver_decl.hpp:332
"Preconditioner" that uses LAPACK's dense LU.
Definition Ifpack2_Details_DenseSolver_decl.hpp:49
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:75
virtual double getComputeTime() const=0
The time (in seconds) spent in compute().
virtual bool isInitialized() const=0
True if the preconditioner has been successfully initialized, else false.
virtual void compute()=0
Set up the numerical values in this preconditioner.
virtual int getNumCompute() const=0
The number of calls to compute().
virtual Teuchos::RCP< const Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > getMatrix() const=0
The input matrix given to the constructor.
virtual int getNumApply() const=0
The number of calls to apply().
virtual double getApplyTime() const=0
The time (in seconds) spent in apply().
virtual void setParameters(const Teuchos::ParameterList &List)=0
Set this preconditioner's parameters.
virtual void apply(const Tpetra::MultiVector< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > &X, Tpetra::MultiVector< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, MatrixType::scalar_type alpha=Teuchos::ScalarTraits< MatrixType::scalar_type >::one(), MatrixType::scalar_type beta=Teuchos::ScalarTraits< MatrixType::scalar_type >::zero()) const=0
Apply the preconditioner to X, putting the result in Y.
virtual Teuchos::RCP< const Tpetra::Map< MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > getRangeMap() const=0
The range Map of this operator.
virtual Teuchos::RCP< const Tpetra::Map< MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > getDomainMap() const=0
The domain Map of this operator.
virtual bool isComputed() const=0
True if the preconditioner has been successfully computed, else false.
virtual void initialize()=0
Set up the graph structure of this preconditioner.
virtual double getInitializeTime() const=0
The time (in seconds) spent in initialize().
virtual int getNumInitialize() const=0
The number of calls to initialize().
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40