43 #ifndef TEUCHOS_SERIALSPDDENSESOLVER_H
44 #define TEUCHOS_SERIALSPDDENSESOLVER_H
59 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
60 #include "Eigen/Dense"
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_);}
342 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_;
409 template<
typename OrdinalType,
typename ScalarType>
413 shouldEquilibrate_(false),
414 equilibratedA_(false),
415 equilibratedB_(false),
418 estimateSolutionErrors_(false),
419 solutionErrorsEstimated_(false),
422 reciprocalConditionEstimated_(false),
423 refineSolution_(false),
424 solutionRefined_(false),
441 template<
typename OrdinalType,
typename ScalarType>
447 template<
typename OrdinalType,
typename ScalarType>
452 reciprocalConditionEstimated_ =
false;
453 solutionRefined_ =
false;
455 solutionErrorsEstimated_ =
false;
456 equilibratedB_ =
false;
460 template<
typename OrdinalType,
typename ScalarType>
464 equilibratedA_ =
false;
482 template<
typename OrdinalType,
typename ScalarType>
488 numRowCols_ =
A->numRows();
497 template<
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!");
520 template<
typename OrdinalType,
typename ScalarType>
523 estimateSolutionErrors_ = flag;
526 refineSolution_ = refineSolution_ || flag;
530 template<
typename OrdinalType,
typename ScalarType>
533 if (factored())
return(0);
536 "SerialSpdDenseSolver<T>::factor: Cannot factor an inverted matrix!");
538 ANORM_ = Matrix_->normOne();
545 if (refineSolution_ ) {
547 AF_ = Factor_->values();
548 LDAF_ = Factor_->stride();
552 if (equilibrate_) ierr = equilibrateMatrix();
554 if (ierr!=0)
return(ierr);
558 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
559 EigenMatrixMap aMap( AF_, numRowCols_, numRowCols_, EigenOuterStride(LDAF_) );
561 if (Matrix_->upper())
567 this->POTRF(Matrix_->UPLO(), numRowCols_, AF_, LDAF_, &INFO_);
577 template<
typename OrdinalType,
typename ScalarType>
591 ierr = equilibrateRHS();
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!");
609 numRowCols_, 1.0, AF_, LDAF_, RHS_->values(), RHS_->stride(), 0.0,
610 LHS_->values(), LHS_->stride());
611 if (INFO_!=0)
return(INFO_);
616 if (!factored()) factor();
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()) );
630 if (Matrix_->upper())
631 lhsMap = lltu_.solve(rhsMap);
633 lhsMap = lltl_.solve(rhsMap);
636 this->POTRS(Matrix_->UPLO(), numRowCols_, RHS_->numCols(), AF_, LDAF_, LHS_->values(), LHS_->stride(), &INFO_);
639 if (INFO_!=0)
return(INFO_);
645 if (shouldEquilibrate() && !equilibratedA_)
646 std::cout <<
"WARNING! SerialSpdDenseSolver<T>::solve: System should be equilibrated!" << std::endl;
649 if (refineSolution_ && !inverted()) ierr1 = applyRefinement();
653 if (equilibrate_) ierr1 = unequilibrateLHS();
658 template<
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_ );
679 this->PORFS(Matrix_->UPLO(), numRowCols_, NRHS, A_, LDA_, AF_, LDAF_,
680 RHS_->values(), RHS_->stride(), LHS_->values(), LHS_->stride(),
681 &FERR_[0], &BERR_[0], &WORK_[0], &PORFS_WORK[0], &INFO_);
683 solutionErrorsEstimated_ =
true;
684 reciprocalConditionEstimated_ =
true;
685 solutionRefined_ =
true;
693 template<
typename OrdinalType,
typename ScalarType>
696 if (R_.size()!=0)
return(0);
698 R_.resize( numRowCols_ );
701 this->POEQU (numRowCols_, AF_, LDAF_, &R_[0], &SCOND_, &AMAX_, &INFO_);
705 shouldEquilibrate_ =
true;
712 template<
typename OrdinalType,
typename ScalarType>
718 if (equilibratedA_)
return(0);
719 if (R_.size()==0) ierr = computeEquilibrateScaling();
720 if (ierr!=0)
return(ierr);
721 if (Matrix_->upper()) {
724 for (j=0; j<numRowCols_; j++) {
726 ScalarType s1 = R_[j];
727 for (i=0; i<=j; i++) {
728 *ptr = *ptr*s1*R_[i];
736 for (j=0; j<numRowCols_; j++) {
738 ptr1 = AF_ + j*LDAF_;
739 ScalarType s1 = R_[j];
740 for (i=0; i<=j; i++) {
741 *ptr = *ptr*s1*R_[i];
743 *ptr1 = *ptr1*s1*R_[i];
752 for (j=0; j<numRowCols_; j++) {
753 ptr = A_ + j + j*LDA_;
754 ScalarType s1 = R_[j];
755 for (i=j; i<numRowCols_; i++) {
756 *ptr = *ptr*s1*R_[i];
764 for (j=0; j<numRowCols_; j++) {
765 ptr = A_ + j + j*LDA_;
766 ptr1 = AF_ + j + j*LDAF_;
767 ScalarType s1 = R_[j];
768 for (i=j; i<numRowCols_; i++) {
769 *ptr = *ptr*s1*R_[i];
771 *ptr1 = *ptr1*s1*R_[i];
778 equilibratedA_ =
true;
785 template<
typename OrdinalType,
typename ScalarType>
791 if (equilibratedB_)
return(0);
792 if (R_.size()==0) ierr = computeEquilibrateScaling();
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++) {
800 for (i=0; i<numRowCols_; i++) {
806 equilibratedB_ =
true;
813 template<
typename OrdinalType,
typename ScalarType>
818 if (!equilibratedB_)
return(0);
820 OrdinalType LDX = LHS_->stride(), NLHS = LHS_->numCols();
821 ScalarType * X = LHS_->values();
823 for (j=0; j<NLHS; j++) {
825 for (i=0; i<numRowCols_; i++) {
836 template<
typename OrdinalType,
typename ScalarType>
839 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
843 if (!factored()) factor();
846 this->POTRI( Matrix_->UPLO(), numRowCols_, AF_, LDAF_, &INFO_);
849 if (Matrix_->upper())
863 template<
typename OrdinalType,
typename ScalarType>
866 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
870 if (reciprocalConditionEstimated()) {
878 if (!factored()) ierr = factor();
879 if (ierr!=0)
return(ierr);
886 std::vector<typename details::lapack_traits<ScalarType>::iwork_type> POCON_WORK( numRowCols_ );
887 this->POCON( Matrix_->UPLO(), numRowCols_, AF_, LDAF_, ANORM_, &RCOND_, &WORK_[0], &POCON_WORK[0], &INFO_);
888 reciprocalConditionEstimated_ =
true;
896 template<
typename OrdinalType,
typename ScalarType>
899 if (Matrix_!=
Teuchos::null) os <<
"Solver Matrix" << std::endl << *Matrix_ << std::endl;
900 if (Factor_!=
Teuchos::null) os <<
"Solver Factored Matrix" << std::endl << *Factor_ << std::endl;
901 if (LHS_ !=
Teuchos::null) os <<
"Solver LHS" << std::endl << *LHS_ << std::endl;
902 if (RHS_ !=
Teuchos::null) os <<
"Solver RHS" << std::endl << *RHS_ << std::endl;