Teuchos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Teuchos_SerialDenseSolver.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Teuchos: Common Tools Package
4 //
5 // Copyright 2004 NTESS and the Teuchos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef _TEUCHOS_SERIALDENSESOLVER_HPP_
11 #define _TEUCHOS_SERIALDENSESOLVER_HPP_
12 
16 #include "Teuchos_CompObject.hpp"
17 #include "Teuchos_BLAS.hpp"
18 #include "Teuchos_LAPACK.hpp"
19 #include "Teuchos_RCP.hpp"
20 #include "Teuchos_ConfigDefs.hpp"
22 #include "Teuchos_ScalarTraits.hpp"
23 
24 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
25 #include "Eigen/Dense"
26 #endif
27 
102 namespace Teuchos {
103 
104  template<typename OrdinalType, typename ScalarType>
105  class SerialDenseSolver : public CompObject, public Object, public BLAS<OrdinalType, ScalarType>,
106  public LAPACK<OrdinalType, ScalarType>
107  {
108 
109  public:
110 
112 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
113  typedef typename Eigen::Matrix<ScalarType,Eigen::Dynamic,Eigen::Dynamic,Eigen::ColMajor> EigenMatrix;
114  typedef typename Eigen::Matrix<ScalarType,Eigen::Dynamic,1> EigenVector;
115  typedef typename Eigen::InnerStride<Eigen::Dynamic> EigenInnerStride;
116  typedef typename Eigen::OuterStride<Eigen::Dynamic> EigenOuterStride;
117  typedef typename Eigen::Map<EigenVector,0,EigenInnerStride > EigenVectorMap;
118  typedef typename Eigen::Map<const EigenVector,0,EigenInnerStride > EigenConstVectorMap;
119  typedef typename Eigen::Map<EigenMatrix,0,EigenOuterStride > EigenMatrixMap;
120  typedef typename Eigen::Map<const EigenMatrix,0,EigenOuterStride > EigenConstMatrixMap;
121  typedef typename Eigen::PermutationMatrix<Eigen::Dynamic,Eigen::Dynamic,OrdinalType> EigenPermutationMatrix;
122  typedef typename Eigen::Array<OrdinalType,Eigen::Dynamic,1> EigenOrdinalArray;
123  typedef typename Eigen::Map<EigenOrdinalArray> EigenOrdinalArrayMap;
124  typedef typename Eigen::Array<ScalarType,Eigen::Dynamic,1> EigenScalarArray;
125  typedef typename Eigen::Map<EigenScalarArray> EigenScalarArrayMap;
126 #endif
127 
129 
132 
134  virtual ~SerialDenseSolver();
136 
138 
139 
142 
144 
150 
152 
153 
155 
157  void factorWithEquilibration(bool flag) {equilibrate_ = flag; return;}
158 
160 
162  void solveWithTranspose(bool flag) {transpose_ = flag; if (flag) TRANS_ = Teuchos::TRANS; else TRANS_ = Teuchos::NO_TRANS; return;}
163 
165 
167  void solveWithTransposeFlag(Teuchos::ETransp trans) {TRANS_ = trans; transpose_ = (trans != Teuchos::NO_TRANS) ? true : false; }
168 
170 
172  void solveToRefinedSolution(bool flag) {refineSolution_ = flag; return;}
173 
175 
179  void estimateSolutionErrors(bool flag);
181 
183 
184 
186 
189  int factor();
190 
192 
195  int solve();
196 
198 
201  int invert();
202 
204 
208 
210 
214  int equilibrateMatrix();
215 
217 
221  int equilibrateRHS();
222 
224 
228  int applyRefinement();
229 
231 
235  int unequilibrateLHS();
236 
238 
246 
248 
249 
251  bool transpose() {return(transpose_);}
252 
254  bool factored() {return(factored_);}
255 
257  bool equilibratedA() {return(equilibratedA_);}
258 
260  bool equilibratedB() {return(equilibratedB_);}
261 
264 
267 
269  bool inverted() {return(inverted_);}
270 
273 
275  bool solved() {return(solved_);}
276 
280 
282 
283 
286 
289 
292 
295 
297  OrdinalType numRows() const {return(M_);}
298 
300  OrdinalType numCols() const {return(N_);}
301 
303  std::vector<OrdinalType> IPIV() const {return(IPIV_);}
304 
306  MagnitudeType ANORM() const {return(ANORM_);}
307 
309  MagnitudeType RCOND() const {return(RCOND_);}
310 
312 
314  MagnitudeType ROWCND() const {return(ROWCND_);}
315 
317 
319  MagnitudeType COLCND() const {return(COLCND_);}
320 
322  MagnitudeType AMAX() const {return(AMAX_);}
323 
325  std::vector<MagnitudeType> FERR() const {return(FERR_);}
326 
328  std::vector<MagnitudeType> BERR() const {return(BERR_);}
329 
331  std::vector<MagnitudeType> R() const {return(R_);}
332 
334  std::vector<MagnitudeType> C() const {return(C_);}
336 
338 
339  void Print(std::ostream& os) const;
342  protected:
343 
344  void allocateWORK() { LWORK_ = 4*N_; WORK_.resize( LWORK_ ); return;}
345  void resetMatrix();
346  void resetVectors();
347 
348 
354  bool factored_;
357  bool solved_;
358  bool inverted_;
362 
364 
365  OrdinalType M_;
366  OrdinalType N_;
367  OrdinalType Min_MN_;
368  OrdinalType LDA_;
369  OrdinalType LDAF_;
370  OrdinalType INFO_;
371  OrdinalType LWORK_;
372 
373  std::vector<OrdinalType> IPIV_;
374 
380 
385 
386  ScalarType * A_;
387  ScalarType * AF_;
388  std::vector<MagnitudeType> FERR_;
389  std::vector<MagnitudeType> BERR_;
390  std::vector<ScalarType> WORK_;
391  std::vector<MagnitudeType> R_;
392  std::vector<MagnitudeType> C_;
393 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
394  Eigen::PartialPivLU<EigenMatrix> lu_;
395 #endif
396 
397  private:
398  // SerialDenseSolver copy constructor (put here because we don't want user access)
399 
402 
403  };
404 
405  namespace details {
406 
407  // Helper traits to distinguish work arrays for real and complex-valued datatypes.
408  template<typename T>
409  struct lapack_traits {
410  typedef int iwork_type;
411  };
412 
413  // Complex-valued specialization
414  template<typename T>
415  struct lapack_traits<std::complex<T> > {
417  };
418 
419  } // end namespace details
420 
421 //=============================================================================
422 
423 template<typename OrdinalType, typename ScalarType>
425  : CompObject(),
426  equilibrate_(false),
427  shouldEquilibrate_(false),
428  equilibratedA_(false),
429  equilibratedB_(false),
430  transpose_(false),
431  factored_(false),
432  estimateSolutionErrors_(false),
433  solutionErrorsEstimated_(false),
434  solved_(false),
435  inverted_(false),
436  reciprocalConditionEstimated_(false),
437  refineSolution_(false),
438  solutionRefined_(false),
439  TRANS_(Teuchos::NO_TRANS),
440  M_(0),
441  N_(0),
442  Min_MN_(0),
443  LDA_(0),
444  LDAF_(0),
445  INFO_(0),
446  LWORK_(0),
447  ANORM_(ScalarTraits<MagnitudeType>::zero()),
448  RCOND_(ScalarTraits<MagnitudeType>::zero()),
449  ROWCND_(ScalarTraits<MagnitudeType>::zero()),
450  COLCND_(ScalarTraits<MagnitudeType>::zero()),
451  AMAX_(ScalarTraits<MagnitudeType>::zero()),
452  A_(0),
453  AF_(0)
454 {
455  resetMatrix();
456 }
457 //=============================================================================
458 
459 template<typename OrdinalType, typename ScalarType>
461 {}
462 
463 //=============================================================================
464 
465 template<typename OrdinalType, typename ScalarType>
467 {
468  LHS_ = Teuchos::null;
469  RHS_ = Teuchos::null;
470  reciprocalConditionEstimated_ = false;
471  solutionRefined_ = false;
472  solved_ = false;
473  solutionErrorsEstimated_ = false;
474  equilibratedB_ = false;
475 }
476 //=============================================================================
477 
478 template<typename OrdinalType, typename ScalarType>
480 {
481  resetVectors();
482  equilibratedA_ = false;
483  factored_ = false;
484  inverted_ = false;
485  M_ = 0;
486  N_ = 0;
487  Min_MN_ = 0;
488  LDA_ = 0;
489  LDAF_ = 0;
495  A_ = 0;
496  AF_ = 0;
497  INFO_ = 0;
498  LWORK_ = 0;
499  R_.resize(0);
500  C_.resize(0);
501 }
502 //=============================================================================
503 
504 template<typename OrdinalType, typename ScalarType>
506 {
507  resetMatrix();
508  Matrix_ = A;
509  Factor_ = A;
510  M_ = A->numRows();
511  N_ = A->numCols();
512  Min_MN_ = TEUCHOS_MIN(M_,N_);
513  LDA_ = A->stride();
514  LDAF_ = LDA_;
515  A_ = A->values();
516  AF_ = A->values();
517  return(0);
518 }
519 //=============================================================================
520 
521 template<typename OrdinalType, typename ScalarType>
524 {
525  // Check that these new vectors are consistent.
526  TEUCHOS_TEST_FOR_EXCEPTION(B->numRows()!=X->numRows() || B->numCols() != X->numCols(), std::invalid_argument,
527  "SerialDenseSolver<T>::setVectors: X and B are not the same size!");
528  TEUCHOS_TEST_FOR_EXCEPTION(B->values()==0, std::invalid_argument,
529  "SerialDenseSolver<T>::setVectors: B is an empty SerialDenseMatrix<T>!");
530  TEUCHOS_TEST_FOR_EXCEPTION(X->values()==0, std::invalid_argument,
531  "SerialDenseSolver<T>::setVectors: X is an empty SerialDenseMatrix<T>!");
532  TEUCHOS_TEST_FOR_EXCEPTION(B->stride()<1, std::invalid_argument,
533  "SerialDenseSolver<T>::setVectors: B has an invalid stride!");
534  TEUCHOS_TEST_FOR_EXCEPTION(X->stride()<1, std::invalid_argument,
535  "SerialDenseSolver<T>::setVectors: X has an invalid stride!");
536 
537  resetVectors();
538  LHS_ = X;
539  RHS_ = B;
540  return(0);
541 }
542 //=============================================================================
543 
544 template<typename OrdinalType, typename ScalarType>
546 {
547  estimateSolutionErrors_ = flag;
548 
549  // If the errors are estimated, this implies that the solution must be refined
550  refineSolution_ = refineSolution_ || flag;
551 }
552 //=============================================================================
553 
554 template<typename OrdinalType, typename ScalarType>
556 
557  if (factored()) return(0); // Already factored
558 
559  TEUCHOS_TEST_FOR_EXCEPTION(inverted(), std::logic_error,
560  "SerialDenseSolver<T>::factor: Cannot factor an inverted matrix!");
561 
562  ANORM_ = Matrix_->normOne(); // Compute 1-Norm of A
563 
564 
565  // If we want to refine the solution, then the factor must
566  // be stored separatedly from the original matrix
567 
568  if (A_ == AF_)
569  if (refineSolution_ ) {
570  Factor_ = rcp( new SerialDenseMatrix<OrdinalType,ScalarType>(Teuchos::Copy, *Matrix_) );
571  AF_ = Factor_->values();
572  LDAF_ = Factor_->stride();
573  }
574 
575  int ierr = 0;
576  if (equilibrate_) ierr = equilibrateMatrix();
577 
578  if (ierr!=0) return(ierr);
579 
580  if ((int)IPIV_.size()!=Min_MN_) IPIV_.resize( Min_MN_ ); // Allocated Pivot vector if not already done.
581 
582  INFO_ = 0;
583 
584 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
585  EigenMatrixMap aMap( AF_, M_, N_, EigenOuterStride(LDAF_) );
586  EigenPermutationMatrix p;
587  EigenOrdinalArray ind;
588  EigenOrdinalArrayMap ipivMap( &IPIV_[0], Min_MN_ );
589 
590  lu_.compute(aMap);
591  p = lu_.permutationP();
592  ind = p.indices();
593 
594  EigenMatrix luMat = lu_.matrixLU();
595  for (OrdinalType j=0; j<aMap.outerSize(); j++) {
596  for (OrdinalType i=0; i<aMap.innerSize(); i++) {
597  aMap(i,j) = luMat(i,j);
598  }
599  }
600  for (OrdinalType i=0; i<ipivMap.innerSize(); i++) {
601  ipivMap(i) = ind(i);
602  }
603 #else
604  this->GETRF(M_, N_, AF_, LDAF_, &IPIV_[0], &INFO_);
605 #endif
606 
607  factored_ = true;
608 
609  return(INFO_);
610 }
611 
612 //=============================================================================
613 
614 template<typename OrdinalType, typename ScalarType>
616 
617  // We will call one of four routines depending on what services the user wants and
618  // whether or not the matrix has been inverted or factored already.
619  //
620  // If the matrix has been inverted, use DGEMM to compute solution.
621  // Otherwise, if the user want the matrix to be equilibrated or wants a refined solution, we will
622  // call the X interface.
623  // Otherwise, if the matrix is already factored we will call the TRS interface.
624  // Otherwise, if the matrix is unfactored we will call the SV interface.
625 
626  int ierr = 0;
627  if (equilibrate_) {
628  ierr = equilibrateRHS();
629  }
630  if (ierr != 0) return(ierr); // Can't equilibrate B, so return.
631 
632  TEUCHOS_TEST_FOR_EXCEPTION( RHS_==Teuchos::null, std::invalid_argument,
633  "SerialDenseSolver<T>::solve: No right-hand side vector (RHS) has been set for the linear system!");
634  TEUCHOS_TEST_FOR_EXCEPTION( LHS_==Teuchos::null, std::invalid_argument,
635  "SerialDenseSolver<T>::solve: No solution vector (LHS) has been set for the linear system!");
636 
637  if (inverted()) {
638 
639  TEUCHOS_TEST_FOR_EXCEPTION( RHS_->values() == LHS_->values(), std::invalid_argument,
640  "SerialDenseSolver<T>::solve: X and B must be different vectors if matrix is inverted.");
641  TEUCHOS_TEST_FOR_EXCEPTION( (equilibratedA_ && !equilibratedB_) || (!equilibratedA_ && equilibratedB_) ,
642  std::logic_error, "SerialDenseSolver<T>::solve: Matrix and vectors must be similarly scaled!");
643 
644  INFO_ = 0;
645  this->GEMM(TRANS_, Teuchos::NO_TRANS, N_, RHS_->numCols(), N_, 1.0, AF_, LDAF_,
646  RHS_->values(), RHS_->stride(), 0.0, LHS_->values(), LHS_->stride());
647  if (INFO_!=0) return(INFO_);
648  solved_ = true;
649  }
650  else {
651 
652  if (!factored()) factor(); // Matrix must be factored
653 
654  TEUCHOS_TEST_FOR_EXCEPTION( (equilibratedA_ && !equilibratedB_) || (!equilibratedA_ && equilibratedB_) ,
655  std::logic_error, "SerialDenseSolver<T>::solve: Matrix and vectors must be similarly scaled!");
656 
657  if (RHS_->values()!=LHS_->values()) {
658  (*LHS_) = (*RHS_); // Copy B to X if needed
659  }
660  INFO_ = 0;
661 
662 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
663  EigenMatrixMap rhsMap( RHS_->values(), RHS_->numRows(), RHS_->numCols(), EigenOuterStride(RHS_->stride()) );
664  EigenMatrixMap lhsMap( LHS_->values(), LHS_->numRows(), LHS_->numCols(), EigenOuterStride(LHS_->stride()) );
665  if ( TRANS_ == Teuchos::NO_TRANS ) {
666  lhsMap = lu_.solve(rhsMap);
667  } else if ( TRANS_ == Teuchos::TRANS ) {
668  lu_.matrixLU().template triangularView<Eigen::Upper>().transpose().solveInPlace(lhsMap);
669  lu_.matrixLU().template triangularView<Eigen::UnitLower>().transpose().solveInPlace(lhsMap);
670  EigenMatrix x = lu_.permutationP().transpose() * lhsMap;
671  lhsMap = x;
672  } else {
673  lu_.matrixLU().template triangularView<Eigen::Upper>().adjoint().solveInPlace(lhsMap);
674  lu_.matrixLU().template triangularView<Eigen::UnitLower>().adjoint().solveInPlace(lhsMap);
675  EigenMatrix x = lu_.permutationP().transpose() * lhsMap;
676  lhsMap = x;
677  }
678 #else
679  this->GETRS(ETranspChar[TRANS_], N_, RHS_->numCols(), AF_, LDAF_, &IPIV_[0], LHS_->values(), LHS_->stride(), &INFO_);
680 #endif
681 
682  if (INFO_!=0) return(INFO_);
683  solved_ = true;
684 
685  }
686 
687  // Warn the user that the matrix should be equilibrated, if it isn't being done already.
688  if (shouldEquilibrate() && !equilibratedA_)
689  std::cout << "WARNING! SerialDenseSolver<T>::solve: System should be equilibrated!" << std::endl;
690 
691  int ierr1=0;
692  if (refineSolution_ && !inverted()) ierr1 = applyRefinement();
693  if (ierr1!=0)
694  return(ierr1);
695 
696  if (equilibrate_) ierr1 = unequilibrateLHS();
697  return(ierr1);
698 }
699 //=============================================================================
700 
701 template<typename OrdinalType, typename ScalarType>
703 {
704  TEUCHOS_TEST_FOR_EXCEPTION(!solved(), std::logic_error,
705  "SerialDenseSolver<T>::applyRefinement: Must have an existing solution!");
706  TEUCHOS_TEST_FOR_EXCEPTION(A_==AF_, std::logic_error,
707  "SerialDenseSolver<T>::applyRefinement: Cannot apply refinement if no original copy of A!");
708 
709 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
710  // Implement templated GERFS or use Eigen.
711  return (-1);
712 #else
713 
714  OrdinalType NRHS = RHS_->numCols();
715  FERR_.resize( NRHS );
716  BERR_.resize( NRHS );
717  allocateWORK();
718 
719  INFO_ = 0;
720  std::vector<typename details::lapack_traits<ScalarType>::iwork_type> GERFS_WORK( N_ );
721  this->GERFS(ETranspChar[TRANS_], N_, NRHS, A_, LDA_, AF_, LDAF_, &IPIV_[0],
722  RHS_->values(), RHS_->stride(), LHS_->values(), LHS_->stride(),
723  &FERR_[0], &BERR_[0], &WORK_[0], &GERFS_WORK[0], &INFO_);
724 
725  solutionErrorsEstimated_ = true;
726  reciprocalConditionEstimated_ = true;
727  solutionRefined_ = true;
728 
729  return(INFO_);
730 #endif
731 }
732 
733 //=============================================================================
734 
735 template<typename OrdinalType, typename ScalarType>
737 {
738  if (R_.size()!=0) return(0); // Already computed
739 
740  R_.resize( M_ );
741  C_.resize( N_ );
742 
743  INFO_ = 0;
744  this->GEEQU (M_, N_, AF_, LDAF_, &R_[0], &C_[0], &ROWCND_, &COLCND_, &AMAX_, &INFO_);
745 
746  if (COLCND_<0.1*ScalarTraits<MagnitudeType>::one() ||
747  ROWCND_<0.1*ScalarTraits<MagnitudeType>::one() ||
749  shouldEquilibrate_ = true;
750 
751  return(INFO_);
752 }
753 
754 //=============================================================================
755 
756 template<typename OrdinalType, typename ScalarType>
758 {
759  OrdinalType i, j;
760  int ierr = 0;
761 
762  if (equilibratedA_) return(0); // Already done.
763  if (R_.size()==0) ierr = computeEquilibrateScaling(); // Compute R and C if needed.
764  if (ierr!=0) return(ierr); // If return value is not zero, then can't equilibrate.
765  if (A_==AF_) {
766  ScalarType * ptr;
767  for (j=0; j<N_; j++) {
768  ptr = A_ + j*LDA_;
769  ScalarType s1 = C_[j];
770  for (i=0; i<M_; i++) {
771  *ptr = *ptr*s1*R_[i];
772  ptr++;
773  }
774  }
775  }
776  else {
777  ScalarType * ptr;
778  ScalarType * ptr1;
779  for (j=0; j<N_; j++) {
780  ptr = A_ + j*LDA_;
781  ptr1 = AF_ + j*LDAF_;
782  ScalarType s1 = C_[j];
783  for (i=0; i<M_; i++) {
784  *ptr = *ptr*s1*R_[i];
785  ptr++;
786  *ptr1 = *ptr1*s1*R_[i];
787  ptr1++;
788  }
789  }
790  }
791 
792  equilibratedA_ = true;
793 
794  return(0);
795 }
796 
797 //=============================================================================
798 
799 template<typename OrdinalType, typename ScalarType>
801 {
802  OrdinalType i, j;
803  int ierr = 0;
804 
805  if (equilibratedB_) return(0); // Already done.
806  if (R_.size()==0) ierr = computeEquilibrateScaling(); // Compute R and C if needed.
807  if (ierr!=0) return(ierr); // Can't count on R and C being computed.
808 
809  MagnitudeType * R_tmp = &R_[0];
810  if (transpose_) R_tmp = &C_[0];
811 
812  OrdinalType LDB = RHS_->stride(), NRHS = RHS_->numCols();
813  ScalarType * B = RHS_->values();
814  ScalarType * ptr;
815  for (j=0; j<NRHS; j++) {
816  ptr = B + j*LDB;
817  for (i=0; i<M_; i++) {
818  *ptr = *ptr*R_tmp[i];
819  ptr++;
820  }
821  }
822 
823  equilibratedB_ = true;
824 
825  return(0);
826 }
827 
828 //=============================================================================
829 
830 template<typename OrdinalType, typename ScalarType>
832 {
833  OrdinalType i, j;
834 
835  if (!equilibratedB_) return(0); // Nothing to do
836 
837  MagnitudeType * C_tmp = &C_[0];
838  if (transpose_) C_tmp = &R_[0];
839 
840  OrdinalType LDX = LHS_->stride(), NLHS = LHS_->numCols();
841  ScalarType * X = LHS_->values();
842  ScalarType * ptr;
843  for (j=0; j<NLHS; j++) {
844  ptr = X + j*LDX;
845  for (i=0; i<N_; i++) {
846  *ptr = *ptr*C_tmp[i];
847  ptr++;
848  }
849  }
850 
851  return(0);
852 }
853 
854 //=============================================================================
855 
856 template<typename OrdinalType, typename ScalarType>
858 {
859 
860  if (!factored()) factor(); // Need matrix factored.
861 
862 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
863  EigenMatrixMap invMap( AF_, M_, N_, EigenOuterStride(LDAF_) );
864  EigenMatrix invMat = lu_.inverse();
865  for (OrdinalType j=0; j<invMap.outerSize(); j++) {
866  for (OrdinalType i=0; i<invMap.innerSize(); i++) {
867  invMap(i,j) = invMat(i,j);
868  }
869  }
870 #else
871 
872  /* This section work with LAPACK Version 3.0 only
873  // Setting LWORK = -1 and calling GETRI will return optimal work space size in WORK_TMP
874  OrdinalType LWORK_TMP = -1;
875  ScalarType WORK_TMP;
876  GETRI( N_, AF_, LDAF_, &IPIV_[0], &WORK_TMP, &LWORK_TMP, &INFO_);
877  LWORK_TMP = WORK_TMP; // Convert to integer
878  if (LWORK_TMP>LWORK_) {
879  LWORK_ = LWORK_TMP;
880  WORK_.resize( LWORK_ );
881  }
882  */
883  // This section will work with any version of LAPACK
884  allocateWORK();
885 
886  INFO_ = 0;
887  this->GETRI( N_, AF_, LDAF_, &IPIV_[0], &WORK_[0], LWORK_, &INFO_);
888 #endif
889 
890  inverted_ = true;
891  factored_ = false;
892 
893  return(INFO_);
894 
895 }
896 
897 //=============================================================================
898 
899 template<typename OrdinalType, typename ScalarType>
901 {
902 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
903  // Implement templated GECON or use Eigen. Eigen currently doesn't have a condition estimation function.
904  return (-1);
905 #else
906  if (reciprocalConditionEstimated()) {
907  Value = RCOND_;
908  return(0); // Already computed, just return it.
909  }
910 
911  if ( ANORM_<ScalarTraits<MagnitudeType>::zero() ) ANORM_ = Matrix_->normOne();
912 
913  int ierr = 0;
914  if (!factored()) ierr = factor(); // Need matrix factored.
915  if (ierr!=0) return(ierr);
916 
917  allocateWORK();
918 
919  // We will assume a one-norm condition number
920  INFO_ = 0;
921  std::vector<typename details::lapack_traits<ScalarType>::iwork_type> GECON_WORK( 2*N_ );
922  this->GECON( '1', N_, AF_, LDAF_, ANORM_, &RCOND_, &WORK_[0], &GECON_WORK[0], &INFO_);
923 
924  reciprocalConditionEstimated_ = true;
925  Value = RCOND_;
926 
927  return(INFO_);
928 #endif
929 }
930 //=============================================================================
931 
932 template<typename OrdinalType, typename ScalarType>
934 
935  if (Matrix_!=Teuchos::null) os << "Solver Matrix" << std::endl << *Matrix_ << std::endl;
936  if (Factor_!=Teuchos::null) os << "Solver Factored Matrix" << std::endl << *Factor_ << std::endl;
937  if (LHS_ !=Teuchos::null) os << "Solver LHS" << std::endl << *LHS_ << std::endl;
938  if (RHS_ !=Teuchos::null) os << "Solver RHS" << std::endl << *RHS_ << std::endl;
939 
940 }
941 
942 } // namespace Teuchos
943 
944 #endif /* _TEUCHOS_SERIALDENSESOLVER_HPP_ */
int equilibrateMatrix()
Equilibrates the this matrix.
void Print(std::ostream &os) const
Print service methods; defines behavior of ostream &lt;&lt; operator.
T magnitudeType
Mandatory typedef for result of magnitude.
Templated serial dense matrix class.
void solveWithTranspose(bool flag)
If flag is true, causes all subsequent function calls to work with the transpose of this matrix...
SerialDenseSolver & operator=(const SerialDenseSolver< OrdinalType, ScalarType > &Source)
std::vector< MagnitudeType > R() const
Returns a pointer to the row scaling vector used for equilibration.
Templated interface class to BLAS routines.
std::vector< MagnitudeType > C_
int solve()
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()..
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Macro for throwing an exception with breakpointing to ease debugging.
bool solved()
Returns true if the current set of vectors has been solved.
bool solutionRefined()
Returns true if the current set of vectors has been refined.
Teuchos header file which uses auto-configuration information to include necessary C++ headers...
Templated BLAS wrapper.
std::vector< MagnitudeType > FERR_
bool reciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getFactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
Object for storing data and providing functionality that is common to all computational classes...
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > Matrix_
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getRHS() const
Returns pointer to current RHS.
std::vector< MagnitudeType > FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
This structure defines some basic traits for a scalar field type.
int unequilibrateLHS()
Unscales the solution vectors if equilibration was used to solve the system.
void solveToRefinedSolution(bool flag)
Causes all solves to compute solution to best ability using iterative refinement. ...
bool factored()
Returns true if matrix is factored (factor available via getFactoredMatrix()).
Functionality and data that is common to all computational classes.
Templated interface class to LAPACK routines.
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getLHS() const
Returns pointer to current LHS.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Deprecated.
bool shouldEquilibrate()
Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system...
bool inverted()
Returns true if matrix inverse has been computed (inverse available via getFactoredMatrix()).
int reciprocalConditionEstimate(MagnitudeType &Value)
Returns the reciprocal of the 1-norm condition number of the this matrix.
The base Teuchos class.
void factorWithEquilibration(bool flag)
Causes equilibration to be called just before the matrix factorization as part of the call to factor...
std::vector< OrdinalType > IPIV_
void solveWithTransposeFlag(Teuchos::ETransp trans)
All subsequent function calls will work with the transpose-type set by this method (Teuchos::NO_TRANS...
int invert()
Inverts the this matrix.
bool transpose()
Returns true if transpose of this matrix has and will be used.
int applyRefinement()
Apply Iterative Refinement.
int computeEquilibrateScaling()
Computes the scaling vector S(i) = 1/sqrt(A(i,i)) of the this matrix.
void estimateSolutionErrors(bool flag)
Causes all solves to estimate the forward and backward solution error.
MagnitudeType ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
bool equilibratedB()
Returns true if RHS is equilibrated (RHS available via getRHS()).
SerialDenseSolver()
Default constructor; matrix should be set using setMatrix(), LHS and RHS set with setVectors()...
bool solutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
The Templated LAPACK Wrapper Class.
std::vector< MagnitudeType > R_
virtual ~SerialDenseSolver()
SerialDenseSolver destructor.
std::vector< OrdinalType > IPIV() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
OrdinalType numRows() const
Returns row dimension of system.
std::vector< MagnitudeType > C() const
Returns a pointer to the column scale vector used for equilibration.
int equilibrateRHS()
Equilibrates the current RHS.
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).
MagnitudeType COLCND() const
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
MagnitudeType ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed)...
std::vector< MagnitudeType > BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
MagnitudeType AMAX() const
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed)...
Defines basic traits for the scalar field type.
TEUCHOSNUMERICS_LIB_DLL_EXPORT const char ETranspChar[]
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
OrdinalType numCols() const
Returns column dimension of system.
Smart reference counting pointer class for automatic garbage collection.
int factor()
Computes the in-place LU factorization of the matrix using the LAPACK routine _GETRF.
std::vector< MagnitudeType > BERR_
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > Factor_
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > RHS_
bool equilibratedA()
Returns true if factor is equilibrated (factor available via getFactoredMatrix()).
#define TEUCHOS_MIN(x, y)
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getMatrix() const
Returns pointer to current matrix.
Reference-counted pointer class and non-member templated function implementations.
static T one()
Returns representation of one for this scalar type.
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > LHS_
A class for solving dense linear problems.
MagnitudeType RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed)...
int setMatrix(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &A)
Sets the pointers for coefficient matrix.
This class creates and provides basic support for dense rectangular matrix of templated type...