43#ifndef TEUCHOS_SERIALSPDDENSESOLVER_H
44#define TEUCHOS_SERIALSPDDENSESOLVER_H
59#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
132 template<
typename OrdinalType,
typename ScalarType>
134 public LAPACK<OrdinalType, ScalarType>
139#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
140 typedef typename Eigen::Matrix<ScalarType,Eigen::Dynamic,Eigen::Dynamic,Eigen::ColMajor> EigenMatrix;
141 typedef typename Eigen::Matrix<ScalarType,Eigen::Dynamic,1> EigenVector;
142 typedef typename Eigen::InnerStride<Eigen::Dynamic> EigenInnerStride;
143 typedef typename Eigen::OuterStride<Eigen::Dynamic> EigenOuterStride;
144 typedef typename Eigen::Map<EigenVector,0,EigenInnerStride > EigenVectorMap;
145 typedef typename Eigen::Map<const EigenVector,0,EigenInnerStride > EigenConstVectorMap;
146 typedef typename Eigen::Map<EigenMatrix,0,EigenOuterStride > EigenMatrixMap;
147 typedef typename Eigen::Map<const EigenMatrix,0,EigenOuterStride > EigenConstMatrixMap;
148 typedef typename Eigen::PermutationMatrix<Eigen::Dynamic,Eigen::Dynamic,OrdinalType> EigenPermutationMatrix;
149 typedef typename Eigen::Array<OrdinalType,Eigen::Dynamic,1> EigenOrdinalArray;
150 typedef typename Eigen::Map<EigenOrdinalArray> EigenOrdinalArrayMap;
151 typedef typename Eigen::Array<ScalarType,Eigen::Dynamic,1> EigenScalarArray;
152 typedef typename Eigen::Map<EigenScalarArray> EigenScalarArrayMap;
330 std::vector<MagnitudeType>
FERR()
const {
return(
FERR_);}
333 std::vector<MagnitudeType>
BERR()
const {
return(
BERR_);}
336 std::vector<MagnitudeType>
R()
const {
return(
R_);}
343 void Print(std::ostream& os)
const;
390 std::vector<MagnitudeType>
R_;
391#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
392 Eigen::LLT<EigenMatrix,Eigen::Upper> lltu_;
393 Eigen::LLT<EigenMatrix,Eigen::Lower> lltl_;
409template<
typename OrdinalType,
typename ScalarType>
441template<
typename OrdinalType,
typename ScalarType>
447template<
typename OrdinalType,
typename ScalarType>
460template<
typename OrdinalType,
typename ScalarType>
482template<
typename OrdinalType,
typename ScalarType>
497template<
typename OrdinalType,
typename ScalarType>
503 "SerialSpdDenseSolver<T>::setVectors: X and B are not the same size!");
505 "SerialSpdDenseSolver<T>::setVectors: B is an empty SerialDenseMatrix<T>!");
507 "SerialSpdDenseSolver<T>::setVectors: X is an empty SerialDenseMatrix<T>!");
509 "SerialSpdDenseSolver<T>::setVectors: B has an invalid stride!");
511 "SerialSpdDenseSolver<T>::setVectors: X has an invalid stride!");
520template<
typename OrdinalType,
typename ScalarType>
530template<
typename OrdinalType,
typename ScalarType>
536 "SerialSpdDenseSolver<T>::factor: Cannot factor an inverted matrix!");
554 if (ierr!=0)
return(ierr);
558#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
577template<
typename OrdinalType,
typename ScalarType>
593 if (ierr != 0)
return(ierr);
596 "SerialSpdDenseSolver<T>::solve: No right-hand side vector (RHS) has been set for the linear system!");
598 "SerialSpdDenseSolver<T>::solve: No solution vector (LHS) has been set for the linear system!");
603 "SerialSpdDenseSolver<T>::solve: X and B must be different vectors if matrix is inverted.");
605 std::logic_error,
"SerialSpdDenseSolver<T>::solve: Matrix and vectors must be similarly scaled!");
619 std::logic_error,
"SerialSpdDenseSolver<T>::solve: Matrix and vectors must be similarly scaled!");
621 if (
RHS_->values()!=
LHS_->values()) {
626#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
627 EigenMatrixMap rhsMap(
RHS_->values(),
RHS_->numRows(),
RHS_->numCols(), EigenOuterStride(
RHS_->stride()) );
628 EigenMatrixMap lhsMap(
LHS_->values(),
LHS_->numRows(),
LHS_->numCols(), EigenOuterStride(
LHS_->stride()) );
631 lhsMap = lltu_.solve(rhsMap);
633 lhsMap = lltl_.solve(rhsMap);
646 std::cout <<
"WARNING! SerialSpdDenseSolver<T>::solve: System should be equilibrated!" << std::endl;
658template<
typename OrdinalType,
typename ScalarType>
662 "SerialSpdDenseSolver<T>::applyRefinement: Must have an existing solution!");
664 "SerialSpdDenseSolver<T>::applyRefinement: Cannot apply refinement if no original copy of A!");
667#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
671 OrdinalType NRHS =
RHS_->numCols();
672 FERR_.resize( NRHS );
673 BERR_.resize( NRHS );
678 std::vector<typename details::lapack_traits<ScalarType>::iwork_type> PORFS_WORK(
numRowCols_ );
693template<
typename OrdinalType,
typename ScalarType>
696 if (
R_.size()!=0)
return(0);
712template<
typename OrdinalType,
typename ScalarType>
720 if (ierr!=0)
return(ierr);
726 ScalarType s1 =
R_[j];
727 for (i=0; i<=j; i++) {
739 ScalarType s1 =
R_[j];
740 for (i=0; i<=j; i++) {
743 *ptr1 = *ptr1*s1*
R_[i];
754 ScalarType s1 =
R_[j];
767 ScalarType s1 =
R_[j];
771 *ptr1 = *ptr1*s1*
R_[i];
785template<
typename OrdinalType,
typename ScalarType>
793 if (ierr!=0)
return(ierr);
795 OrdinalType LDB =
RHS_->stride(), NRHS =
RHS_->numCols();
796 ScalarType *
B =
RHS_->values();
798 for (j=0; j<NRHS; j++) {
813template<
typename OrdinalType,
typename ScalarType>
820 OrdinalType LDX =
LHS_->stride(), NLHS =
LHS_->numCols();
821 ScalarType * X =
LHS_->values();
823 for (j=0; j<NLHS; j++) {
836template<
typename OrdinalType,
typename ScalarType>
839#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
863template<
typename OrdinalType,
typename ScalarType>
866#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
879 if (ierr!=0)
return(ierr);
886 std::vector<typename details::lapack_traits<ScalarType>::iwork_type> POCON_WORK(
numRowCols_ );
896template<
typename OrdinalType,
typename ScalarType>
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.
Templated serial dense matrix class.
Templated class for solving dense linear problems.
Templated serial, dense, symmetric matrix class.
BLAS(void)
Default constructor.
CompObject()
Default constructor.
void GEMM(ETransp transa, ETransp transb, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
General matrix-matrix multiply.
void PORFS(const char &UPLO, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const ScalarType *AF, const OrdinalType &ldaf, 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 system of linear equations when the coefficient matrix is symmetr...
void POCON(const char &UPLO, const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, const ScalarType &anorm, ScalarType *rcond, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Estimates the reciprocal of the condition number (1-norm) of a real symmetric positive definite matri...
void POTRS(const char &UPLO, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Solves a system of linear equations A*X=B, where A is a symmetric positive definite matrix factored b...
void POEQU(const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, MagnitudeType *S, MagnitudeType *scond, MagnitudeType *amax, OrdinalType *info) const
Computes row and column scalings intended to equilibrate a symmetric positive definite matrix A and r...
void POTRI(const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *info) const
Computes the inverse of a real symmetric positive definite matrix A using the Cholesky factorization ...
LAPACK(void)
Default Constructor.
void POTRF(const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *info) const
Computes Cholesky factorization of a real symmetric positive definite matrix A.
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 dense rectangular matrix of templated type.
SerialSpdDenseSolver(const SerialSpdDenseSolver< OrdinalType, ScalarType > &Source)
MagnitudeType SCOND()
Ratio of smallest to largest equilibration scale factors for the this matrix (returns -1 if not yet c...
int setMatrix(const RCP< SerialSymDenseMatrix< OrdinalType, ScalarType > > &A_in)
Sets the pointers for coefficient matrix.
std::vector< MagnitudeType > BERR_
std::vector< MagnitudeType > FERR_
OrdinalType numCols() const
Returns column dimension of system.
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
RCP< SerialSymDenseMatrix< OrdinalType, ScalarType > > Factor_
SerialSpdDenseSolver & operator=(const SerialSpdDenseSolver< OrdinalType, ScalarType > &Source)
std::vector< MagnitudeType > FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
int applyRefinement()
Apply Iterative Refinement.
bool estimateSolutionErrors_
bool transpose()
Returns true if transpose of this matrix has and will be used.
bool reciprocalConditionEstimated_
int solve()
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()....
std::vector< MagnitudeType > R() const
Returns a pointer to the row scaling vector used for equilibration.
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).
int factor()
Computes the in-place Cholesky factorization of the matrix using the LAPACK routine DPOTRF.
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getRHS() const
Returns pointer to current RHS.
void estimateSolutionErrors(bool flag)
Causes all solves to estimate the forward and backward solution error.
bool solutionRefined()
Returns true if the current set of vectors has been refined.
virtual ~SerialSpdDenseSolver()
SerialSpdDenseSolver destructor.
RCP< SerialSymDenseMatrix< OrdinalType, ScalarType > > getMatrix() const
Returns pointer to current matrix.
bool solved()
Returns true if the current set of vectors has been solved.
bool reciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
RCP< SerialSymDenseMatrix< OrdinalType, ScalarType > > getFactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
void solveToRefinedSolution(bool flag)
Causes all solves to compute solution to best ability using iterative refinement.
OrdinalType numRows() const
Returns row dimension of system.
std::vector< ScalarType > WORK_
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > RHS_
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > LHS_
std::vector< int > IWORK_
MagnitudeType ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
bool solutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
MagnitudeType AMAX() const
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed).
bool factored()
Returns true if matrix is factored (factor available via AF() and LDAF()).
void Print(std::ostream &os) const
Print service methods; defines behavior of ostream << operator.
SerialSpdDenseSolver()
Default constructor; matrix should be set using setMatrix(), LHS and RHS set with setVectors().
std::vector< MagnitudeType > BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
bool solutionErrorsEstimated_
int unequilibrateLHS()
Unscales the solution vectors if equilibration was used to solve the system.
void factorWithEquilibration(bool flag)
Causes equilibration to be called just before the matrix factorization as part of the call to factor.
int equilibrateMatrix()
Equilibrates the this matrix.
int invert()
Inverts the this matrix.
int reciprocalConditionEstimate(MagnitudeType &Value)
Returns the reciprocal of the 1-norm condition number of the this matrix.
bool equilibratedB()
Returns true if RHS is equilibrated (RHS available via B() and LDB()).
RCP< SerialSymDenseMatrix< OrdinalType, ScalarType > > Matrix_
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getLHS() const
Returns pointer to current LHS.
bool inverted()
Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()).
std::vector< MagnitudeType > R_
int computeEquilibrateScaling()
Computes the scaling vector S(i) = 1/sqrt(A(i,i) of the this matrix.
MagnitudeType RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed).
int equilibrateRHS()
Equilibrates the current RHS.
bool shouldEquilibrate()
Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system.
bool equilibratedA()
Returns true if factor is equilibrated (factor available via AF() and LDAF()).
This class creates and provides basic support for symmetric, positive-definite dense matrices of temp...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Macro for throwing an exception with breakpointing to ease debugging.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Deprecated.
Teuchos implementation details.
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).