42#ifndef _TEUCHOS_SERIALQRDENSESOLVER_HPP_
43#define _TEUCHOS_SERIALQRDENSESOLVER_HPP_
57#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
131 template<
typename OrdinalType,
typename ScalarType>
133 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;
320 std::vector<ScalarType>
tau()
const {
return(
TAU_);}
330 void Print(std::ostream& os)
const;
382#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
383 Eigen::HouseholderQR<EigenMatrix> qr_;
399template<
typename OrdinalType,
typename ScalarType>
435template<
typename OrdinalType,
typename ScalarType>
441template<
typename OrdinalType,
typename ScalarType>
452template<
typename OrdinalType,
typename ScalarType>
481template<
typename OrdinalType,
typename ScalarType>
485 "SerialQRDenseSolver<T>::setMatrix: the matrix A must have A.numRows() >= A.numCols()!");
508template<
typename OrdinalType,
typename ScalarType>
514 "SerialQRDenseSolver<T>::setVectors: X and B have different numbers of columns!");
516 "SerialQRDenseSolver<T>::setVectors: B is an empty SerialDenseMatrix<T>!");
518 "SerialQRDenseSolver<T>::setVectors: X is an empty SerialDenseMatrix<T>!");
520 "SerialQRDenseSolver<T>::setVectors: B has an invalid stride!");
522 "SerialQRDenseSolver<T>::setVectors: X has an invalid stride!");
532template<
typename OrdinalType,
typename ScalarType>
540 if (ierr!=0)
return(ierr);
548#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
549 EigenMatrixMap aMap(
AF_,
M_,
N_, EigenOuterStride(
LDAF_) );
550 EigenScalarArray
tau;
554 for (OrdinalType i=0; i<tauMap.innerSize(); i++) {
557 EigenMatrix qrMat = qr_.matrixQR();
558 for (OrdinalType j=0; j<aMap.outerSize(); j++) {
559 for (OrdinalType i=0; i<aMap.innerSize(); i++) {
560 aMap(i,j) = qrMat(i,j);
574template<
typename OrdinalType,
typename ScalarType>
588 if (ierr != 0)
return(ierr);
591 "SerialQRDenseSolver<T>::solve: No right-hand side vector (RHS) has been set for the linear system!");
593 "SerialQRDenseSolver<T>::solve: No solution vector (LHS) has been set for the linear system!");
596 "SerialQRDenseSolver<T>::solve: No transpose: must have LHS_->numRows() = N_");
599 "SerialQRDenseSolver<T>::solve: Transpose: must have LHS_->numRows() = M_");
603 std::cout <<
"WARNING! SerialQRDenseSolver<T>::solve: System should be equilibrated!" << std::endl;
609 std::logic_error,
"SerialQRDenseSolver<T>::solve: Matrix and vectors must be similarly scaled!");
612 for (OrdinalType j=0; j<
RHS_->numCols(); j++) {
613 for (OrdinalType i=0; i<
RHS_->numRows(); i++) {
614 (*TMP_)(i,j) = (*
RHS_)(i,j);
631 for (OrdinalType j = 0; j <
RHS_->numCols(); j++ ) {
632 for (OrdinalType i =
N_; i <
M_; i++ ) {
639 for (OrdinalType j = 0; j <
LHS_->numCols(); j++ ) {
640 for (OrdinalType i = 0; i <
LHS_->numRows(); i++ ) {
641 (*LHS_)(i, j) = (*
TMP_)(i,j);
660template<
typename OrdinalType,
typename ScalarType>
675 for (j = 0; j <
N_; j++) {
676 for (i = 0; i <
M_; i++) {
686 }
else if (
ANORM_ > bignum) {
689 }
else if (
ANORM_ == mZero) {
699template<
typename OrdinalType,
typename ScalarType>
717 for (j = 0; j <
N_; j++) {
718 for (i = 0; i <
M_; i++) {
729 }
else if (
ANORM_ > bignum) {
733 }
else if (
ANORM_ == mZero) {
745template<
typename OrdinalType,
typename ScalarType>
763 for (j = 0; j <
RHS_->numCols(); j++) {
764 for (i = 0; i <
RHS_->numRows(); i++) {
777 }
else if (
BNORM_ > bignum) {
791template<
typename OrdinalType,
typename ScalarType>
828template<
typename OrdinalType,
typename ScalarType>
844#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
845 EigenMatrixMap qMap(
Q_,
M_,
N_, EigenOuterStride(
LDQ_) );
846 EigenMatrix qMat = qr_.householderQ();
847 for (OrdinalType j=0; j<qMap.outerSize(); j++) {
848 for (OrdinalType i=0; i<qMap.innerSize(); i++) {
849 qMap(i,j) = qMat(i,j);
865template<
typename OrdinalType,
typename ScalarType>
879 for (OrdinalType j=0; j<
N_; j++) {
880 for (OrdinalType i=0; i<=j; i++) {
881 (*FactorR_)(i,j) = (*
Factor_)(i,j);
892template<
typename OrdinalType,
typename ScalarType>
898 "SerialQRDenseSolver<T>::multiplyQ: C.numRows() != M_!");
900 "SerialQRDenseSolver<T>::multiplyQ: C is an empty SerialDenseMatrix<T>!");
908#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
909 EigenMatrixMap cMap(
C.values(),
C.numRows(),
C.numCols(), EigenOuterStride(
C.stride()) );
912 cMap = qr_.householderQ() * cMap;
915 cMap = qr_.householderQ().adjoint() * cMap;
916 for (OrdinalType j = 0; j <
C.numCols(); j++) {
917 for (OrdinalType i =
N_; i <
C.numRows(); i++ ) {
939 for (OrdinalType j = 0; j <
C.numCols(); j++) {
940 for (OrdinalType i =
N_; i <
C.numRows(); i++ ) {
953template<
typename OrdinalType,
typename ScalarType>
959 "SerialQRDenseSolver<T>::solveR: must have N_ < C.numRows() < M_!");
961 "SerialQRDenseSolver<T>::solveR: C is an empty SerialDenseMatrix<T>!");
969#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
970 EigenMatrixMap cMap(
C.values(),
N_,
C.numCols(), EigenOuterStride(
C.stride()) );
972 for (OrdinalType j=0; j<
N_; j++) {
980 qr_.matrixQR().topRows(
N_).template triangularView<Eigen::Upper>().solveInPlace(cMap);
983 qr_.matrixQR().topRows(
N_).template triangularView<Eigen::Upper>().adjoint().solveInPlace(cMap);
998 RMAT, LDRMAT,
C.values(),
C.stride(), &
INFO_);
1004 RMAT, LDRMAT,
C.values(),
C.stride(), &
INFO_);
1015template<
typename OrdinalType,
typename ScalarType>
1020 if (
Q_!=
Teuchos::null) os <<
"Solver Factor Q" << std::endl << *
Q_ << 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.
#define TEUCHOS_MIN(x, y)
#define TEUCHOS_MAX(x, y)
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 class for solving dense linear problems.
BLAS(void)
Default constructor.
CompObject()
Default constructor.
void TRTRS(const char &UPLO, const char &TRANS, const char &DIAG, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Solves a triangular linear system of the form A*X=B or A**T*X=B, where A is a triangular matrix.
void UNGQR(const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Compute explicit QR factor from QR factorization (GEQRF) (complex case).
void UNMQR(const char &SIDE, const char &TRANS, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, const ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *C, const OrdinalType &ldc, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Apply Householder reflectors (complex case).
void LASCL(const char &TYPE, const OrdinalType &kl, const OrdinalType &ku, const MagnitudeType cfrom, const MagnitudeType cto, const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *info) const
Multiplies the m by n matrix A by the real scalar cto/cfrom.
LAPACK(void)
Default Constructor.
void GEQRF(const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Computes a QR factorization of a general m by n matrix A.
Object(int tracebackModeIn=-1)
Default Constructor.
Smart reference counting pointer class for automatic garbage collection.
This class creates and provides basic support for dense rectangular matrix of templated type.
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > FactorR_
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > LHS_
void factorWithEquilibration(bool flag)
Causes equilibration to be called just before the matrix factorization as part of the call to factor.
SerialQRDenseSolver & operator=(const SerialQRDenseSolver< OrdinalType, ScalarType > &Source)
int solveR(ETransp transr, SerialDenseMatrix< OrdinalType, ScalarType > &C)
Solve input matrix on the left with the upper triangular matrix R or its adjoint.
bool solved()
Returns true if the current set of vectors has been solved.
std::vector< ScalarType > WORK_
bool equilibratedB()
Returns true if RHS is equilibrated (RHS available via getRHS()).
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getLHS() const
Returns pointer to current LHS.
bool formedR()
Returns true if R has been formed explicitly.
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > FactorQ_
int formR()
Explicitly forms the upper triangular matrix R.
ScalarTraits< ScalarType >::magnitudeType MagnitudeType
bool shouldEquilibrate()
Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the 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).
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > Matrix_
SerialQRDenseSolver(const SerialQRDenseSolver< OrdinalType, ScalarType > &Source)
OrdinalType equilibrationModeA_
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getMatrix() const
Returns pointer to current matrix.
void solveWithTranspose(bool flag)
If flag is true, causes all subsequent function calls to work with the adjoint of this matrix,...
bool transpose()
Returns true if adjoint of this matrix has and will be used.
int multiplyQ(ETransp transq, SerialDenseMatrix< OrdinalType, ScalarType > &C)
Left multiply the input matrix by the unitary matrix Q or its adjoint.
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getQ() const
Returns pointer to Q (assuming factorization has been performed).
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > TMP_
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > RHS_
MagnitudeType ANORM() const
Returns the absolute value of the largest element of this matrix (returns -1 if not yet computed).
OrdinalType numCols() const
Returns column dimension of system.
bool formedQ()
Returns true if Q has been formed explicitly.
std::vector< ScalarType > tau() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
OrdinalType equilibrationModeB_
OrdinalType numRows() const
Returns row dimension of system.
virtual ~SerialQRDenseSolver()
SerialQRDenseSolver destructor.
int factor()
Computes the in-place QR factorization of the matrix using the LAPACK routine _GETRF or the Eigen cla...
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > Factor_
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getFactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
int solve()
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()....
SerialQRDenseSolver()
Default constructor; matrix should be set using setMatrix(), LHS and RHS set with setVectors().
void solveWithTransposeFlag(Teuchos::ETransp trans)
All subsequent function calls will work with the transpose-type set by this method (Teuchos::NO_TRANS...
bool factored()
Returns true if matrix is factored (factor available via getFactoredMatrix()).
void Print(std::ostream &os) const
Print service methods; defines behavior of ostream << operator.
std::vector< ScalarType > TAU_
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getR() const
Returns pointer to R (assuming factorization has been performed).
int setMatrix(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &A)
Sets the pointers for coefficient matrix.
int formQ()
Explicitly forms the unitary matrix Q.
int equilibrateMatrix()
Equilibrates the this matrix.
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getRHS() const
Returns pointer to current RHS.
int equilibrateRHS()
Equilibrates the current RHS.
bool equilibratedA()
Returns true if factor is equilibrated (factor available via getFactoredMatrix()).
int computeEquilibrateScaling()
Determines if this matrix should be scaled.
int unequilibrateLHS()
Unscales the solution vectors if equilibration was used to solve the system.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Macro for throwing an exception with breakpointing to ease debugging.
TEUCHOSNUMERICS_LIB_DLL_EXPORT const char ETranspChar[]
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Deprecated.
TEUCHOSNUMERICS_LIB_DLL_EXPORT const char EUploChar[]
TEUCHOSNUMERICS_LIB_DLL_EXPORT const char ETypeChar[]
TEUCHOSNUMERICS_LIB_DLL_EXPORT const char EDiagChar[]
TEUCHOSNUMERICS_LIB_DLL_EXPORT const char ESideChar[]
Teuchos implementation details.
This structure defines some basic traits for a scalar field type.
static magnitudeType magnitude(T a)
Returns the magnitudeType of the scalar type a.
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 const bool isComplex
Determines if scalar type is std::complex.
static magnitudeType sfmin()
Returns safe minimum (sfmin), such that 1/sfmin does not overflow.
static magnitudeType prec()
Returns eps*base.