42 #ifndef _TEUCHOS_SERIALQRDENSESOLVER_HPP_
43 #define _TEUCHOS_SERIALQRDENSESOLVER_HPP_
57 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
58 #include "Eigen/Dense"
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_);}
323 MagnitudeType
ANORM()
const {
return(ANORM_);}
329 void Print(std::ostream& os)
const;
334 void allocateWORK() { LWORK_ = 4*N_; WORK_.resize( LWORK_ );
return;}
340 bool shouldEquilibrate_;
343 OrdinalType equilibrationModeA_;
344 OrdinalType equilibrationModeB_;
364 std::vector<ScalarType> TAU_;
366 MagnitudeType ANORM_;
367 MagnitudeType BNORM_;
369 RCP<SerialDenseMatrix<OrdinalType, ScalarType> > Matrix_;
370 RCP<SerialDenseMatrix<OrdinalType, ScalarType> > LHS_;
371 RCP<SerialDenseMatrix<OrdinalType, ScalarType> > TMP_;
372 RCP<SerialDenseMatrix<OrdinalType, ScalarType> > RHS_;
373 RCP<SerialDenseMatrix<OrdinalType, ScalarType> > Factor_;
374 RCP<SerialDenseMatrix<OrdinalType, ScalarType> > FactorQ_;
375 RCP<SerialDenseMatrix<OrdinalType, ScalarType> > FactorR_;
381 std::vector<ScalarType> WORK_;
382 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
383 Eigen::HouseholderQR<EigenMatrix> qr_;
390 SerialQRDenseSolver & operator=(
const SerialQRDenseSolver<OrdinalType, ScalarType>& Source);
399 template<
typename OrdinalType,
typename ScalarType>
403 shouldEquilibrate_(false),
404 equilibratedA_(false),
405 equilibratedB_(false),
406 equilibrationModeA_(0),
407 equilibrationModeB_(0),
411 matrixIsZero_(false),
435 template<
typename OrdinalType,
typename ScalarType>
441 template<
typename OrdinalType,
typename ScalarType>
444 LHS_ = Teuchos::null;
445 TMP_ = Teuchos::null;
446 RHS_ = Teuchos::null;
448 equilibratedB_ =
false;
452 template<
typename OrdinalType,
typename ScalarType>
453 void SerialQRDenseSolver<OrdinalType,ScalarType>::resetMatrix()
456 equilibratedA_ =
false;
457 equilibrationModeA_ = 0;
458 equilibrationModeB_ = 0;
460 matrixIsZero_ =
false;
481 template<
typename OrdinalType,
typename ScalarType>
485 "SerialQRDenseSolver<T>::setMatrix: the matrix A must have A.numRows() >= A.numCols()!");
494 Min_MN_ = TEUCHOS_MIN(M_,N_);
508 template<
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!");
532 template<
typename OrdinalType,
typename ScalarType>
535 if (factored())
return(0);
539 if (equilibrate_) ierr = equilibrateMatrix();
540 if (ierr!=0)
return(ierr);
543 if ((
int)TAU_.size()!=Min_MN_) TAU_.resize( Min_MN_ );
548 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
549 EigenMatrixMap aMap( AF_, M_, N_, EigenOuterStride(LDAF_) );
550 EigenScalarArray tau;
551 EigenScalarArrayMap tauMap(&TAU_[0],TEUCHOS_MIN(M_,N_));
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);
564 this->GEQRF(M_, N_, AF_, LDAF_, &TAU_[0], &WORK_[0], LWORK_, &INFO_);
574 template<
typename OrdinalType,
typename ScalarType>
586 ierr = equilibrateRHS();
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_");
602 if (shouldEquilibrate() && !equilibratedA_)
603 std::cout <<
"WARNING! SerialQRDenseSolver<T>::solve: System should be equilibrated!" << std::endl;
606 if (!factored()) factor();
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);
649 ierr = unequilibrateLHS();
660 template<
typename OrdinalType,
typename ScalarType>
675 for (j = 0; j < N_; j++) {
676 for (i = 0; i < M_; i++) {
683 if (ANORM_ > mZero && ANORM_ < smlnum) {
685 shouldEquilibrate_ =
true;
686 }
else if (ANORM_ > bignum) {
688 shouldEquilibrate_ =
true;
689 }
else if (ANORM_ == mZero) {
691 matrixIsZero_ =
true;
699 template<
typename OrdinalType,
typename ScalarType>
702 if (equilibratedA_)
return(0);
717 for (j = 0; j < N_; j++) {
718 for (i = 0; i < M_; i++) {
725 if (ANORM_ > mZero && ANORM_ < smlnum) {
727 this->LASCL(Teuchos::ETypeChar[
Teuchos::FULL], BW, BW, ANORM_, smlnum, M_, N_, A_, LDA_, &INFO_);
728 equilibrationModeA_ = 1;
729 }
else if (ANORM_ > bignum) {
731 this->LASCL(Teuchos::ETypeChar[
Teuchos::FULL], BW, BW, ANORM_, bignum, M_, N_, A_, LDA_, &INFO_);
732 equilibrationModeA_ = 2;
733 }
else if (ANORM_ == mZero) {
735 matrixIsZero_ =
true;
738 equilibratedA_ =
true;
745 template<
typename OrdinalType,
typename ScalarType>
748 if (equilibratedB_)
return(0);
763 for (j = 0; j <RHS_->numCols(); j++) {
764 for (i = 0; i < RHS_->numRows(); i++) {
772 if (BNORM_ > mZero && BNORM_ < smlnum) {
774 this->LASCL(Teuchos::ETypeChar[
Teuchos::FULL], BW, BW, BNORM_, smlnum, RHS_->numRows(), RHS_->numCols(),
775 RHS_->values(), RHS_->stride(), &INFO_);
776 equilibrationModeB_ = 1;
777 }
else if (BNORM_ > bignum) {
779 this->LASCL(Teuchos::ETypeChar[
Teuchos::FULL], BW, BW, BNORM_, bignum, RHS_->numRows(), RHS_->numCols(),
780 RHS_->values(), RHS_->stride(), &INFO_);
781 equilibrationModeB_ = 2;
784 equilibratedB_ =
true;
791 template<
typename OrdinalType,
typename ScalarType>
794 if (!equilibratedB_)
return(0);
808 if (equilibrationModeA_ == 1) {
809 this->LASCL(Teuchos::ETypeChar[
Teuchos::FULL], BW, BW, ANORM_, smlnum, LHS_->numRows(), LHS_->numCols(),
810 LHS_->values(), LHS_->stride(), &INFO_);
811 }
else if (equilibrationModeA_ == 2) {
812 this->LASCL(Teuchos::ETypeChar[
Teuchos::FULL], BW, BW, ANORM_, bignum, LHS_->numRows(), LHS_->numCols(),
813 LHS_->values(), LHS_->stride(), &INFO_);
815 if (equilibrationModeB_ == 1) {
816 this->LASCL(Teuchos::ETypeChar[
Teuchos::FULL], BW, BW, smlnum, BNORM_, LHS_->numRows(), LHS_->numCols(),
817 LHS_->values(), LHS_->stride(), &INFO_);
818 }
else if (equilibrationModeB_ == 2) {
819 this->LASCL(Teuchos::ETypeChar[
Teuchos::FULL], BW, BW, bignum, BNORM_, LHS_->numRows(), LHS_->numCols(),
820 LHS_->values(), LHS_->stride(), &INFO_);
828 template<
typename OrdinalType,
typename ScalarType>
832 if (!factored()) factor();
837 Q_ = FactorQ_->values();
838 LDQ_ = FactorQ_->stride();
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);
853 this->UNGQR(M_, N_, N_, Q_, LDQ_, &TAU_[0], &WORK_[0], LWORK_, &INFO_);
856 if (INFO_!=0)
return(INFO_);
865 template<
typename OrdinalType,
typename ScalarType>
869 if (!factored()) factor();
874 R_ = FactorR_->values();
875 LDR_ = FactorR_->stride();
879 for (OrdinalType j=0; j<N_; j++) {
880 for (OrdinalType i=0; i<=j; i++) {
881 (*FactorR_)(i,j) = (*Factor_)(i,j);
892 template<
typename OrdinalType,
typename ScalarType>
898 "SerialQRDenseSolver<T>::multiplyQ: C.numRows() != M_!");
900 "SerialQRDenseSolver<T>::multiplyQ: C is an empty SerialDenseMatrix<T>!");
903 if (!factored()) factor();
908 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
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++ ) {
930 this->UNMQR(ESideChar[SIDE], ETranspChar[
NO_TRANS], M_, C.
numCols(), N_, AF_, LDAF_,
931 &TAU_[0], C.
values(), C.
stride(), &WORK_[0], LWORK_, &INFO_);
936 this->UNMQR(ESideChar[SIDE], ETranspChar[
TRANS], M_, C.
numCols(), N_, AF_, LDAF_,
937 &TAU_[0], C.
values(), C.
stride(), &WORK_[0], LWORK_, &INFO_);
939 for (OrdinalType j = 0; j < C.
numCols(); j++) {
940 for (OrdinalType i = N_; i < C.
numRows(); i++ ) {
953 template<
typename OrdinalType,
typename ScalarType>
959 "SerialQRDenseSolver<T>::solveR: must have N_ < C.numRows() < M_!");
961 "SerialQRDenseSolver<T>::solveR: C is an empty SerialDenseMatrix<T>!");
964 if (!factored()) factor();
969 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
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);
991 ScalarType* RMAT = (formedR_) ? R_ : AF_;
992 OrdinalType LDRMAT = (formedR_) ? LDR_ : LDAF_;
997 this->TRTRS(EUploChar[UPLO], ETranspChar[
NO_TRANS], EDiagChar[DIAG], N_, C.
numCols(),
1003 this->TRTRS(EUploChar[UPLO], ETranspChar[
TRANS], EDiagChar[DIAG], N_, C.
numCols(),
1015 template<
typename OrdinalType,
typename ScalarType>
1018 if (Matrix_!=Teuchos::null) os <<
"Solver Matrix" << std::endl << *Matrix_ << std::endl;
1019 if (Factor_!=Teuchos::null) os <<
"Solver Factored Matrix" << std::endl << *Factor_ << std::endl;
1020 if (Q_!=Teuchos::null) os <<
"Solver Factor Q" << std::endl << *Q_ << std::endl;
1021 if (LHS_ !=Teuchos::null) os <<
"Solver LHS" << std::endl << *LHS_ << std::endl;
1022 if (RHS_ !=Teuchos::null) os <<
"Solver RHS" << std::endl << *RHS_ << std::endl;