10#ifndef _TEUCHOS_SERIALBANDDENSESOLVER_HPP_
11#define _TEUCHOS_SERIALBANDDENSESOLVER_HPP_
27#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
133 template<
typename OrdinalType,
typename ScalarType>
135 public LAPACK<OrdinalType, ScalarType>
141#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
142 typedef typename Eigen::Matrix<ScalarType,Eigen::Dynamic,Eigen::Dynamic,Eigen::ColMajor> EigenMatrix;
143 typedef typename Eigen::internal::BandMatrix<ScalarType,Eigen::Dynamic,Eigen::Dynamic,Eigen::ColMajor> EigenBandMatrix;
144 typedef typename Eigen::Matrix<ScalarType,Eigen::Dynamic,1> EigenVector;
145 typedef typename Eigen::InnerStride<Eigen::Dynamic> EigenInnerStride;
146 typedef typename Eigen::OuterStride<Eigen::Dynamic> EigenOuterStride;
147 typedef typename Eigen::Map<EigenVector,0,EigenInnerStride > EigenVectorMap;
148 typedef typename Eigen::Map<const EigenVector,0,EigenInnerStride > EigenConstVectorMap;
149 typedef typename Eigen::Map<EigenMatrix,0,EigenOuterStride > EigenMatrixMap;
150 typedef typename Eigen::Map<const EigenMatrix,0,EigenOuterStride > EigenConstMatrixMap;
151 typedef typename Eigen::PermutationMatrix<Eigen::Dynamic,Eigen::Dynamic,OrdinalType> EigenPermutationMatrix;
152 typedef typename Eigen::Array<OrdinalType,Eigen::Dynamic,1> EigenOrdinalArray;
153 typedef typename Eigen::Map<EigenOrdinalArray> EigenOrdinalArrayMap;
154 typedef typename Eigen::Array<ScalarType,Eigen::Dynamic,1> EigenScalarArray;
155 typedef typename Eigen::Map<EigenScalarArray> EigenScalarArrayMap;
322 std::vector<OrdinalType>
IPIV()
const {
return(IPIV_);}
325 MagnitudeType
ANORM()
const {
return(ANORM_);}
328 MagnitudeType
RCOND()
const {
return(RCOND_);}
333 MagnitudeType
ROWCND()
const {
return(ROWCND_);}
338 MagnitudeType
COLCND()
const {
return(COLCND_);}
341 MagnitudeType
AMAX()
const {
return(AMAX_);}
344 std::vector<MagnitudeType>
FERR()
const {
return(FERR_);}
347 std::vector<MagnitudeType>
BERR()
const {
return(BERR_);}
350 std::vector<MagnitudeType>
R()
const {
return(R_);}
353 std::vector<MagnitudeType>
C()
const {
return(C_);}
359 void Print(std::ostream& os)
const;
363 void allocateWORK() { LWORK_ = 3*N_; WORK_.resize( LWORK_ );
return;}
369 bool shouldEquilibrate_;
374 bool estimateSolutionErrors_;
375 bool solutionErrorsEstimated_;
377 bool reciprocalConditionEstimated_;
378 bool refineSolution_;
379 bool solutionRefined_;
393 std::vector<OrdinalType> IPIV_;
394 std::vector<int> IWORK_;
396 MagnitudeType ANORM_;
397 MagnitudeType RCOND_;
398 MagnitudeType ROWCND_;
399 MagnitudeType COLCND_;
402 RCP<SerialBandDenseMatrix<OrdinalType, ScalarType> > Matrix_;
403 RCP<SerialDenseMatrix<OrdinalType, ScalarType> > LHS_;
404 RCP<SerialDenseMatrix<OrdinalType, ScalarType> > RHS_;
405 RCP<SerialBandDenseMatrix<OrdinalType, ScalarType> > Factor_;
409 std::vector<MagnitudeType> FERR_;
410 std::vector<MagnitudeType> BERR_;
411 std::vector<ScalarType> WORK_;
412 std::vector<MagnitudeType> R_;
413 std::vector<MagnitudeType> C_;
414#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
415 Eigen::PartialPivLU<EigenMatrix> lu_;
427 using namespace details;
431template<
typename OrdinalType,
typename ScalarType>
435 shouldEquilibrate_(false),
436 equilibratedA_(false),
437 equilibratedB_(false),
440 estimateSolutionErrors_(false),
441 solutionErrorsEstimated_(false),
443 reciprocalConditionEstimated_(false),
444 refineSolution_(false),
445 solutionRefined_(false),
468template<
typename OrdinalType,
typename ScalarType>
474template<
typename OrdinalType,
typename ScalarType>
475void SerialBandDenseSolver<OrdinalType,ScalarType>::resetVectors()
477 LHS_ = Teuchos::null;
478 RHS_ = Teuchos::null;
479 reciprocalConditionEstimated_ =
false;
480 solutionRefined_ =
false;
482 solutionErrorsEstimated_ =
false;
483 equilibratedB_ =
false;
487template<
typename OrdinalType,
typename ScalarType>
488void SerialBandDenseSolver<OrdinalType,ScalarType>::resetMatrix()
491 equilibratedA_ =
false;
513template<
typename OrdinalType,
typename ScalarType>
517 OrdinalType m = AB->numRows();
518 OrdinalType n = AB->numCols();
519 OrdinalType kl = AB->lowerBandwidth();
520 OrdinalType ku = AB->upperBandwidth();
524 "SerialBandDenseSolver<T>::setMatrix: A is an empty SerialBandDenseMatrix<T>!");
533 Min_MN_ = TEUCHOS_MIN(M_,N_);
543template<
typename OrdinalType,
typename ScalarType>
549 "SerialBandDenseSolver<T>::setVectors: X and B are not the same size!");
551 "SerialBandDenseSolver<T>::setVectors: B is an empty SerialDenseMatrix<T>!");
553 "SerialBandDenseSolver<T>::setVectors: X is an empty SerialDenseMatrix<T>!");
555 "SerialBandDenseSolver<T>::setVectors: B has an invalid stride!");
557 "SerialBandDenseSolver<T>::setVectors: X has an invalid stride!");
566template<
typename OrdinalType,
typename ScalarType>
569 estimateSolutionErrors_ = flag;
572 refineSolution_ = refineSolution_ || flag;
576template<
typename OrdinalType,
typename ScalarType>
581 ANORM_ = Matrix_->normOne();
587 if (refineSolution_ ) {
589 AF_ = Factor_->values();
590 LDAF_ = Factor_->stride();
596 if (ierr!=0)
return(ierr);
598 if (IPIV_.size()==0) IPIV_.resize( N_ );
602#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
603 EigenMatrixMap aMap( AF_+KL_, KL_+KU_+1, N_, EigenOuterStride(LDAF_) );
604 EigenMatrix fullMatrix(M_, N_);
605 for (OrdinalType j=0; j<N_; j++) {
606 for (OrdinalType i=TEUCHOS_MAX(0,j-KU_); i<=TEUCHOS_MIN(M_-1, j+KL_); i++) {
607 fullMatrix(i,j) = aMap(KU_+i-j, j);
611 EigenPermutationMatrix p;
612 EigenOrdinalArray ind;
613 EigenOrdinalArrayMap ipivMap( &IPIV_[0], Min_MN_ );
615 lu_.compute(fullMatrix);
616 p = lu_.permutationP();
619 for (OrdinalType i=0; i<ipivMap.innerSize(); i++) {
624 this->
GBTRF(M_, N_, KL_, KU_, AF_, LDAF_, &IPIV_[0], &INFO_);
634template<
typename OrdinalType,
typename ScalarType>
646 if (ierr != 0)
return(ierr);
649 "SerialBandDenseSolver<T>::solve: No right-hand side vector (RHS) has been set for the linear system!");
651 "SerialBandDenseSolver<T>::solve: No solution vector (LHS) has been set for the linear system!");
656 std::logic_error,
"SerialBandDenseSolver<T>::solve: Matrix and vectors must be similarly scaled!");
659 std::cout <<
"WARNING! SerialBandDenseSolver<T>::solve: System should be equilibrated!" << std::endl;
661 if (RHS_->values()!=LHS_->values()) {
666#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
667 EigenMatrixMap rhsMap( RHS_->values(), RHS_->numRows(), RHS_->numCols(), EigenOuterStride(RHS_->stride()) );
668 EigenMatrixMap lhsMap( LHS_->values(), LHS_->numRows(), LHS_->numCols(), EigenOuterStride(LHS_->stride()) );
670 lhsMap = lu_.solve(rhsMap);
672 lu_.matrixLU().template triangularView<Eigen::Upper>().transpose().solveInPlace(lhsMap);
673 lu_.matrixLU().template triangularView<Eigen::UnitLower>().transpose().solveInPlace(lhsMap);
674 EigenMatrix x = lu_.permutationP().transpose() * lhsMap;
677 lu_.matrixLU().template triangularView<Eigen::Upper>().adjoint().solveInPlace(lhsMap);
678 lu_.matrixLU().template triangularView<Eigen::UnitLower>().adjoint().solveInPlace(lhsMap);
679 EigenMatrix x = lu_.permutationP().transpose() * lhsMap;
683 this->
GBTRS(ETranspChar[TRANS_], N_, KL_, KU_, RHS_->numCols(), AF_, LDAF_, &IPIV_[0], LHS_->values(), LHS_->stride(), &INFO_);
686 if (INFO_!=0)
return(INFO_);
700template<
typename OrdinalType,
typename ScalarType>
704 "SerialBandDenseSolver<T>::applyRefinement: Must have an existing solution!");
706 "SerialBandDenseSolver<T>::applyRefinement: Cannot apply refinement if no original copy of A!");
708#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
713 OrdinalType NRHS = RHS_->numCols();
714 FERR_.resize( NRHS );
715 BERR_.resize( NRHS );
719 std::vector<typename details::lapack_traits<ScalarType>::iwork_type> GBRFS_WORK( N_ );
720 this->
GBRFS(ETranspChar[TRANS_], N_, KL_, KU_, NRHS, A_+KL_, LDA_, AF_, LDAF_, &IPIV_[0],
721 RHS_->values(), RHS_->stride(), LHS_->values(), LHS_->stride(),
722 &FERR_[0], &BERR_[0], &WORK_[0], &GBRFS_WORK[0], &INFO_);
724 solutionErrorsEstimated_ =
true;
725 reciprocalConditionEstimated_ =
true;
726 solutionRefined_ =
true;
734template<
typename OrdinalType,
typename ScalarType>
737 if (R_.size()!=0)
return(0);
743 this->
GBEQU (M_, N_, KL_, KU_, AF_+KL_, LDAF_, &R_[0], &C_[0], &ROWCND_, &COLCND_, &AMAX_, &INFO_);
748 shouldEquilibrate_ =
true;
755template<
typename OrdinalType,
typename ScalarType>
761 if (equilibratedA_)
return(0);
763 if (ierr!=0)
return(ierr);
768 for (j=0; j<N_; j++) {
769 ptr = A_ + KL_ + j*LDA_ + TEUCHOS_MAX(KU_-j, 0);
770 ScalarType s1 = C_[j];
771 for (i=TEUCHOS_MAX(0,j-KU_); i<=TEUCHOS_MIN(M_-1,j+KL_); i++) {
781 for (j=0; j<N_; j++) {
782 ptr = A_ + KL_ + j*LDA_ + TEUCHOS_MAX(KU_-j, 0);
783 ScalarType s1 = C_[j];
784 for (i=TEUCHOS_MAX(0,j-KU_); i<=TEUCHOS_MIN(M_-1,j+KL_); i++) {
789 for (j=0; j<N_; j++) {
790 ptrU = AF_ + j*LDAF_ + TEUCHOS_MAX(KL_+KU_-j, 0);
791 ptrL = AF_ + KL_ + KU_ + 1 + j*LDAF_;
792 ScalarType s1 = C_[j];
793 for (i=TEUCHOS_MAX(0,j-(KL_+KU_)); i<=TEUCHOS_MIN(M_-1,j); i++) {
794 *ptrU = *ptrU*s1*R_[i];
797 for (i=TEUCHOS_MAX(0,j); i<=TEUCHOS_MIN(M_-1,j+KL_); i++) {
798 *ptrL = *ptrL*s1*R_[i];
804 equilibratedA_ =
true;
811template<
typename OrdinalType,
typename ScalarType>
817 if (equilibratedB_)
return(0);
819 if (ierr!=0)
return(ierr);
821 MagnitudeType * R_tmp = &R_[0];
822 if (transpose_) R_tmp = &C_[0];
824 OrdinalType LDB = RHS_->stride(), NRHS = RHS_->numCols();
825 ScalarType * B = RHS_->values();
827 for (j=0; j<NRHS; j++) {
829 for (i=0; i<M_; i++) {
835 equilibratedB_ =
true;
843template<
typename OrdinalType,
typename ScalarType>
848 if (!equilibratedB_)
return(0);
850 MagnitudeType * C_tmp = &C_[0];
851 if (transpose_) C_tmp = &R_[0];
853 OrdinalType LDX = LHS_->stride(), NLHS = LHS_->numCols();
854 ScalarType * X = LHS_->values();
856 for (j=0; j<NLHS; j++) {
858 for (i=0; i<N_; i++) {
869template<
typename OrdinalType,
typename ScalarType>
872#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
886 if (ierr!=0)
return(ierr);
892 std::vector<typename details::lapack_traits<ScalarType>::iwork_type> GBCON_WORK( N_ );
893 this->
GBCON(
'1', N_, KL_, KU_, AF_, LDAF_, &IPIV_[0], ANORM_, &RCOND_, &WORK_[0], &GBCON_WORK[0], &INFO_);
894 reciprocalConditionEstimated_ =
true;
902template<
typename OrdinalType,
typename ScalarType>
905 if (Matrix_!=Teuchos::null) os <<
"Solver Matrix" << std::endl << *Matrix_ << std::endl;
906 if (Factor_!=Teuchos::null) os <<
"Solver Factored Matrix" << std::endl << *Factor_ << std::endl;
907 if (LHS_ !=Teuchos::null) os <<
"Solver LHS" << std::endl << *LHS_ << std::endl;
908 if (RHS_ !=Teuchos::null) os <<
"Solver RHS" << std::endl << *RHS_ << std::endl;
Templated interface class to BLAS routines.
Object for storing data and providing functionality that is common to all computational classes.
Teuchos header file which uses auto-configuration information to include necessary C++ headers.
Templated interface class to LAPACK routines.
Reference-counted pointer class and non-member templated function implementations.
Defines basic traits for the scalar field type.
Templated serial dense matrix class.
Templated serial dense matrix class.
Templated class for solving dense linear problems.
BLAS(void)
Default constructor.
CompObject()
Default constructor.
void GBRFS(const char &TRANS, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const ScalarType *AF, const OrdinalType &ldaf, const OrdinalType *IPIV, const ScalarType *B, const OrdinalType &ldb, ScalarType *X, const OrdinalType &ldx, ScalarType *FERR, ScalarType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Improves the computed solution to a banded system of linear equations and provides error bounds and b...
void GBTRF(const OrdinalType &m, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, OrdinalType *info) const
Computes an LU factorization of a general banded m by n matrix A using partial pivoting with row inte...
void GBEQU(const OrdinalType &m, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const ScalarType *A, const OrdinalType &lda, MagnitudeType *R, MagnitudeType *C, MagnitudeType *rowcond, MagnitudeType *colcond, MagnitudeType *amax, OrdinalType *info) const
Computes row and column scalings intended to equilibrate an m by n banded matrix A and reduce its con...
void GBCON(const char &NORM, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const ScalarType *A, const OrdinalType &lda, const OrdinalType *IPIV, const ScalarType &anorm, ScalarType *rcond, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Estimates the reciprocal of the condition number of a general banded real matrix A,...
void GBTRS(const char &TRANS, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Solves a system of linear equations A*X=B or A'*X=B with a general banded n by n matrix A using the L...
LAPACK(void)
Default Constructor.
Object(int tracebackModeIn=-1)
Default Constructor.
Ptr< T > ptr(T *p)
Create a pointer to an object from a raw pointer.
Smart reference counting pointer class for automatic garbage collection.
This class creates and provides basic support for banded dense matrices of templated type.
bool equilibratedA()
Returns true if factor is equilibrated (factor available via AF() and LDAF()).
void estimateSolutionErrors(bool flag)
Causes all solves to estimate the forward and backward solution error.
void solveWithTranspose(bool flag)
MagnitudeType ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
bool shouldEquilibrate()
Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system.
MagnitudeType ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed).
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getRHS() const
Returns pointer to current RHS.
OrdinalType numRows() const
Returns row dimension of system.
int setVectors(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &X, const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &B)
Sets the pointers for left and right hand side vector(s).
bool factored()
Returns true if matrix is factored (factor available via AF() and LDAF()).
std::vector< OrdinalType > IPIV() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
int applyRefinement()
Apply Iterative Refinement.
bool solutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
OrdinalType numLower() const
Returns lower bandwidth of system.
bool solved()
Returns true if the current set of vectors has been solved.
MagnitudeType RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed).
OrdinalType numUpper() const
Returns upper bandwidth of system.
virtual ~SerialBandDenseSolver()
SerialBandDenseSolver destructor.
int setMatrix(const RCP< SerialBandDenseMatrix< OrdinalType, ScalarType > > &AB)
Sets the pointer for coefficient matrix.
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getLHS() const
Returns pointer to current LHS.
std::vector< MagnitudeType > FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
void factorWithEquilibration(bool flag)
int factor()
Computes the in-place LU factorization of the matrix using the LAPACK routine _GBTRF.
std::vector< MagnitudeType > BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
std::vector< MagnitudeType > R() const
Returns a pointer to the row scaling vector used for equilibration.
bool transpose()
Returns true if transpose of this matrix has and will be used.
MagnitudeType AMAX() const
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed).
std::vector< MagnitudeType > C() const
Returns a pointer to the column scale vector used for equilibration.
MagnitudeType COLCND() const
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
int computeEquilibrateScaling()
Computes the scaling vector S(i) = 1/sqrt(A(i,i)) of the this matrix.
RCP< SerialBandDenseMatrix< OrdinalType, ScalarType > > getFactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
bool reciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
int reciprocalConditionEstimate(MagnitudeType &Value)
Returns the reciprocal of the 1-norm condition number of the this matrix.
int unequilibrateLHS()
Unscales the solution vectors if equilibration was used to solve the system.
bool equilibratedB()
Returns true if RHS is equilibrated (RHS available via B() and LDB()).
void Print(std::ostream &os) const
Print service methods; defines behavior of ostream << operator.
int equilibrateMatrix()
Equilibrates the this matrix.
bool solutionRefined()
Returns true if the current set of vectors has been refined.
OrdinalType numCols() const
Returns column dimension of system.
int equilibrateRHS()
Equilibrates the current RHS.
void solveWithTransposeFlag(Teuchos::ETransp trans)
RCP< SerialBandDenseMatrix< OrdinalType, ScalarType > > getMatrix() const
Returns pointer to current matrix.
int solve()
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()....
void solveToRefinedSolution(bool flag)
Set whether or not to use iterative refinement to improve solutions to linear systems.
SerialBandDenseSolver()
Default constructor; matrix should be set using setMatrix(), LHS and RHS set with setVectors().
This class creates and provides basic support for dense rectangular matrix of templated type.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Macro for throwing an exception with breakpointing to ease debugging.
The Teuchos namespace contains all of the classes, structs and enums used by Teuchos,...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Deprecated.
This structure defines some basic traits for a scalar field type.
static magnitudeType rmax()
Overflow theshold - (base^emax)*(1-eps).
static T one()
Returns representation of one for this scalar type.
T magnitudeType
Mandatory typedef for result of magnitude.
static T zero()
Returns representation of zero for this scalar type.
static magnitudeType rmin()
Returns the underflow threshold - base^(emin-1).