10#ifndef IFPACK2_DETAILS_TRIDISOLVER_DEF_HPP
11#define IFPACK2_DETAILS_TRIDISOLVER_DEF_HPP
13#include "Ifpack2_LocalFilter.hpp"
14#include "Teuchos_LAPACK.hpp"
15#include "Tpetra_MultiVector.hpp"
16#include "Tpetra_Map.hpp"
17#include "Tpetra_Import.hpp"
18#include "Tpetra_Export.hpp"
22#include "Teuchos_DefaultMpiComm.hpp"
24#include "Teuchos_DefaultSerialComm.hpp"
34template <
class MatrixType>
36 TriDiSolver(
const Teuchos::RCP<const row_matrix_type>& A)
38 , initializeTime_(0.0)
44 , isInitialized_(
false)
45 , isComputed_(
false) {}
47template <
class MatrixType>
48Teuchos::RCP<const typename TriDiSolver<MatrixType, false>::map_type>
51 A_.is_null(), std::runtime_error,
52 "Ifpack2::Details::TriDiSolver::"
53 "getDomainMap: The input matrix A is null. Please call setMatrix() with a "
54 "nonnull input matrix before calling this method.");
57 return A_->getRangeMap();
60template <
class MatrixType>
61Teuchos::RCP<const typename TriDiSolver<MatrixType, false>::map_type>
64 A_.is_null(), std::runtime_error,
65 "Ifpack2::Details::TriDiSolver::"
66 "getRangeMap: The input matrix A is null. Please call setMatrix() with a "
67 "nonnull input matrix before calling this method.");
70 return A_->getDomainMap();
73template <
class MatrixType>
79template <
class MatrixType>
81 return isInitialized_;
84template <
class MatrixType>
89template <
class MatrixType>
91 return numInitialize_;
94template <
class MatrixType>
99template <
class MatrixType>
104template <
class MatrixType>
107 return initializeTime_;
110template <
class MatrixType>
116template <
class MatrixType>
122template <
class MatrixType>
123Teuchos::RCP<const typename TriDiSolver<MatrixType, false>::row_matrix_type>
128template <
class MatrixType>
131 isInitialized_ =
false;
133 A_local_ = Teuchos::null;
134 A_local_tridi_.reshape(0);
138template <
class MatrixType>
140 setMatrix(
const Teuchos::RCP<const row_matrix_type>& A) {
144 const global_size_t numRows = A->getRangeMap()->getGlobalNumElements();
145 const global_size_t
numCols = A->getDomainMap()->getGlobalNumElements();
147 numRows !=
numCols, std::invalid_argument,
148 "Ifpack2::Details::TriDiSolver::"
149 "setMatrix: Input matrix must be (globally) square. "
150 "The matrix you provided is "
151 << numRows <<
" by " <<
numCols <<
".");
160template <
class MatrixType>
167 using Teuchos::TimeMonitor;
168 const std::string
timerName(
"Ifpack2::Details::TriDiSolver::initialize");
171 if (
timer.is_null()) {
181 A_.is_null(), std::runtime_error,
182 "Ifpack2::Details::TriDiSolver::"
183 "initialize: The input matrix A is null. Please call setMatrix() "
184 "with a nonnull input before calling this method.");
187 !A_->hasColMap(), std::invalid_argument,
188 "Ifpack2::Details::TriDiSolver: "
189 "The constructor's input matrix must have a column Map, "
190 "so that it has local indices.");
196 if (A_->getComm()->getSize() > 1) {
203 A_local_.is_null(), std::logic_error,
204 "Ifpack2::Details::TriDiSolver::"
205 "initialize: A_local_ is null after it was supposed to have been "
206 "initialized. Please report this bug to the Ifpack2 developers.");
209 const size_t numRows = A_local_->getLocalNumRows();
210 const size_t numCols = A_local_->getLocalNumCols();
212 numRows !=
numCols, std::logic_error,
213 "Ifpack2::Details::TriDiSolver::"
214 "initialize: Local filter matrix is not square. This should never happen. "
215 "Please report this bug to the Ifpack2 developers.");
216 A_local_tridi_.reshape(numRows);
217 ipiv_.resize(numRows);
218 std::fill(ipiv_.begin(), ipiv_.end(), 0);
220 isInitialized_ =
true;
227template <
class MatrixType>
230template <
class MatrixType>
233 const std::string
timerName(
"Ifpack2::Details::TriDiSolver::compute");
236 if (
timer.is_null()) {
246 A_.is_null(), std::runtime_error,
247 "Ifpack2::Details::TriDiSolver::"
248 "compute: The input matrix A is null. Please call setMatrix() with a "
249 "nonnull input, then call initialize(), before calling this method.");
252 A_local_.is_null(), std::logic_error,
253 "Ifpack2::Details::TriDiSolver::"
254 "compute: A_local_ is null. Please report this bug to the Ifpack2 "
258 if (!this->isInitialized()) {
261 extract(A_local_tridi_, *A_local_);
263 factor(A_local_tridi_, ipiv_());
271template <
class MatrixType>
273 const Teuchos::ArrayView<int>& ipiv) {
275 std::fill(ipiv.begin(), ipiv.end(), 0);
277 Teuchos::LAPACK<int, scalar_type>
lapack;
279 lapack.GTTRF(A.numRowsCols(),
284 ipiv.getRawPtr(), &
INFO);
287 INFO < 0, std::logic_error,
288 "Ifpack2::Details::TriDiSolver::factor: "
289 "LAPACK's _GTTRF (tridiagonal LU factorization with partial pivoting) "
290 "was called incorrectly. INFO = "
292 "Please report this bug to the Ifpack2 developers.");
297 INFO > 0, std::runtime_error,
298 "Ifpack2::Details::TriDiSolver::factor: "
299 "LAPACK's _GTTRF (tridiagonal LU factorization with partial pivoting) "
300 "reports that the computed U factor is exactly singular. U("
301 <<
INFO <<
"," <<
INFO <<
") (one-based index i) is exactly zero. This probably "
302 "means that the input matrix has a singular diagonal block.");
305template <
class MatrixType>
306void TriDiSolver<MatrixType, false>::
307 applyImpl(
const MV& X,
309 const Teuchos::ETransp mode,
310 const scalar_type alpha,
311 const scalar_type beta)
const {
312 using Teuchos::ArrayRCP;
313 using Teuchos::CONJ_TRANS;
316 using Teuchos::rcpFromRef;
317 using Teuchos::TRANS;
319 const int numVecs =
static_cast<int>(X.getNumVectors());
320 if (alpha == STS::zero()) {
321 if (beta == STS::zero()) {
325 Y.putScalar(STS::zero());
327 Y.scale(STS::zero());
330 Teuchos::LAPACK<int, scalar_type> lapack;
336 if (beta == STS::zero() && Y.isConstantStride() && alpha == STS::one()) {
338 Y_tmp = rcpFromRef(Y);
340 Y_tmp = rcp(
new MV(createCopy(X)));
341 if (alpha != STS::one()) {
345 const int Y_stride =
static_cast<int>(Y_tmp->getStride());
346 ArrayRCP<scalar_type> Y_view = Y_tmp->get1dViewNonConst();
347 scalar_type*
const Y_ptr = Y_view.getRawPtr();
350 (mode == CONJ_TRANS ?
'C' : (mode ==
TRANS ?
'T' :
'N'));
351 lapack.GTTRS(trans, A_local_tridi_.numRowsCols(), numVecs,
355 A_local_tridi_.DU2(),
356 ipiv_.getRawPtr(), Y_ptr, Y_stride, &INFO);
357 TEUCHOS_TEST_FOR_EXCEPTION(
358 INFO != 0, std::runtime_error,
359 "Ifpack2::Details::TriDiSolver::"
360 "applyImpl: LAPACK's _GTTRS (tridiagonal solve using LU factorization "
361 "with partial pivoting) failed with INFO = "
362 << INFO <<
" != 0.");
364 if (beta != STS::zero()) {
365 Y.update(alpha, *Y_tmp, beta);
366 }
else if (!Y.isConstantStride()) {
367 deep_copy(Y, *Y_tmp);
372template <
class MatrixType>
374 apply(
const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
375 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
376 Teuchos::ETransp
mode,
379 using Teuchos::ArrayView;
383 using Teuchos::rcpFromRef;
385 const std::string
timerName(
"Ifpack2::Details::TriDiSolver::apply");
387 if (
timer.is_null()) {
398 !isComputed_, std::runtime_error,
399 "Ifpack2::Details::TriDiSolver::apply: "
400 "You must have called the compute() method before you may call apply(). "
401 "You may call the apply() method as many times as you want after calling "
402 "compute() once, but you must have called compute() at least once.");
404 const size_t numVecs = X.getNumVectors();
407 numVecs != Y.getNumVectors(), std::runtime_error,
408 "Ifpack2::Details::TriDiSolver::apply: X and Y have different numbers "
410 << X.getNumVectors() <<
", but Y has "
411 << X.getNumVectors() <<
".");
420 const bool multipleProcs = (A_->getRowMap()->getComm()->getSize() >= 1);
425 X_local = X.offsetView(A_local_->getDomainMap(), 0);
426 Y_local = Y.offsetViewNonConst(A_local_->getRangeMap(), 0);
443template <
class MatrixType>
445 std::ostringstream
out;
450 out <<
"\"Ifpack2::Details::TriDiSolver\": ";
455 out <<
"Initialized: " << (isInitialized() ?
"true" :
"false") <<
", "
456 <<
"Computed: " << (isComputed() ?
"true" :
"false") <<
", ";
459 out <<
"Matrix: null";
461 out <<
"Matrix: not null"
462 <<
", Global matrix dimensions: ["
463 << A_->getGlobalNumRows() <<
", " << A_->getGlobalNumCols() <<
"]";
470template <
class MatrixType>
472 const Teuchos::EVerbosityLevel
verbLevel)
const {
474 using Teuchos::FancyOStream;
475 using Teuchos::OSTab;
477 using Teuchos::rcpFromRef;
482 RCP<FancyOStream> ptrOut = rcpFromRef(out);
484 if (this->getObjectLabel() !=
"") {
485 out <<
"label: " << this->getObjectLabel() << endl;
487 out <<
"initialized: " << (isInitialized_ ?
"true" :
"false") << endl
488 <<
"computed: " << (isComputed_ ?
"true" :
"false") << endl
489 <<
"number of initialize calls: " << numInitialize_ << endl
490 <<
"number of compute calls: " << numCompute_ << endl
491 <<
"number of apply calls: " << numApply_ << endl
492 <<
"total time in seconds in initialize: " << initializeTime_ << endl
493 <<
"total time in seconds in compute: " << computeTime_ << endl
494 <<
"total time in seconds in apply: " << applyTime_ << endl;
495 if (verbLevel >= Teuchos::VERB_EXTREME) {
496 out <<
"A_local_tridi_:" << endl;
497 A_local_tridi_.print(out);
499 out <<
"ipiv_: " << Teuchos::toString(ipiv_) << endl;
503template <
class MatrixType>
505 const Teuchos::EVerbosityLevel
verbLevel)
const {
507 using Teuchos::FancyOStream;
508 using Teuchos::OSTab;
510 using Teuchos::rcpFromRef;
520 out <<
"Ifpack2::Details::TriDiSolver:" <<
endl;
526 const Teuchos::Comm<int>&
comm = *(A_->getRowMap()->getComm());
530 out <<
"Ifpack2::Details::TriDiSolver:" <<
endl;
545template <
class MatrixType>
547 const row_matrix_type&
A_local) {
548 using Teuchos::Array;
549 using Teuchos::ArrayView;
550 typedef local_ordinal_type LO;
551 typedef typename Teuchos::ArrayView<LO>::size_type size_type;
567 static_cast<size_type
>(
A_local.getLocalMaxNumRowEntries());
592 const scalar_type val = values[
k];
606template <
class MatrixType>
611template <
class MatrixType>
612Teuchos::RCP<const typename TriDiSolver<MatrixType, true>::map_type>
617template <
class MatrixType>
618Teuchos::RCP<const typename TriDiSolver<MatrixType, true>::map_type>
623template <
class MatrixType>
628template <
class MatrixType>
633template <
class MatrixType>
638template <
class MatrixType>
643template <
class MatrixType>
648template <
class MatrixType>
653template <
class MatrixType>
659template <
class MatrixType>
665template <
class MatrixType>
671template <
class MatrixType>
672Teuchos::RCP<const typename TriDiSolver<MatrixType, true>::row_matrix_type>
677template <
class MatrixType>
682template <
class MatrixType>
687template <
class MatrixType>
692template <
class MatrixType>
697template <
class MatrixType>
699 const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
700 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
701 Teuchos::ETransp
mode,
707template <
class MatrixType>
713template <
class MatrixType>
715 const Teuchos::EVerbosityLevel
verbLevel)
const {
722#define IFPACK2_DETAILS_TRIDISOLVER_INSTANT(S, LO, GO, N) \
723 template class Ifpack2::Details::TriDiSolver<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.
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:75
MatrixType::scalar_type scalar_type
The type of entries in the input (global) matrix.
Definition Ifpack2_Details_TriDiSolver_decl.hpp:72
MatrixType::scalar_type scalar_type
The type of entries in the input (global) matrix.
Definition Ifpack2_Details_TriDiSolver_decl.hpp:331
"Preconditioner" that uses LAPACK's tridi LU.
Definition Ifpack2_Details_TriDiSolver_decl.hpp:50
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