Teuchos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Teuchos_SerialBandDenseSolver.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_SERIALBANDDENSESOLVER_HPP_
11 #define _TEUCHOS_SERIALBANDDENSESOLVER_HPP_
12 
17 #include "Teuchos_CompObject.hpp"
18 #include "Teuchos_BLAS.hpp"
19 #include "Teuchos_LAPACK.hpp"
20 #include "Teuchos_RCP.hpp"
21 #include "Teuchos_ConfigDefs.hpp"
25 #include "Teuchos_ScalarTraits.hpp"
26 
27 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
28 #include "Eigen/Dense"
29 #endif
30 
31 namespace Teuchos {
32 
133  template<typename OrdinalType, typename ScalarType>
134  class SerialBandDenseSolver : public CompObject, public Object, public BLAS<OrdinalType, ScalarType>,
135  public LAPACK<OrdinalType, ScalarType>
136  {
137 
138  public:
139 
141 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
142  typedef typename Eigen::Matrix<ScalarType,Eigen::Dynamic,Eigen::Dynamic,Eigen::ColMajor> EigenMatrix;
143  typedef typename Eigen::internal::BandMatrix<ScalarType,Eigen::Dynamic,Eigen::Dynamic,Eigen::ColMajor> EigenBandMatrix;
144  typedef typename Eigen::Matrix<ScalarType,Eigen::Dynamic,1> EigenVector;
145  typedef typename Eigen::InnerStride<Eigen::Dynamic> EigenInnerStride;
146  typedef typename Eigen::OuterStride<Eigen::Dynamic> EigenOuterStride;
147  typedef typename Eigen::Map<EigenVector,0,EigenInnerStride > EigenVectorMap;
148  typedef typename Eigen::Map<const EigenVector,0,EigenInnerStride > EigenConstVectorMap;
149  typedef typename Eigen::Map<EigenMatrix,0,EigenOuterStride > EigenMatrixMap;
150  typedef typename Eigen::Map<const EigenMatrix,0,EigenOuterStride > EigenConstMatrixMap;
151  typedef typename Eigen::PermutationMatrix<Eigen::Dynamic,Eigen::Dynamic,OrdinalType> EigenPermutationMatrix;
152  typedef typename Eigen::Array<OrdinalType,Eigen::Dynamic,1> EigenOrdinalArray;
153  typedef typename Eigen::Map<EigenOrdinalArray> EigenOrdinalArrayMap;
154  typedef typename Eigen::Array<ScalarType,Eigen::Dynamic,1> EigenScalarArray;
155  typedef typename Eigen::Map<EigenScalarArray> EigenScalarArrayMap;
156 #endif
157 
159 
160 
163 
165  virtual ~SerialBandDenseSolver();
166 
168 
170 
173 
175 
180 
182 
184 
187  void factorWithEquilibration(bool flag) {equilibrate_ = flag; return;}
188 
192  void solveWithTranspose(bool flag) {transpose_ = flag; if (flag) TRANS_ = Teuchos::TRANS; else TRANS_ = Teuchos::NO_TRANS; return;}
193 
196  void solveWithTransposeFlag(Teuchos::ETransp trans) {TRANS_ = trans; transpose_ = (trans != Teuchos::NO_TRANS) ? true : false;}
197 
199  void solveToRefinedSolution(bool flag) {refineSolution_ = flag; return;}
200 
202 
205  void estimateSolutionErrors(bool flag);
207 
209 
210 
212 
215  int factor();
216 
218 
221  int solve();
222 
224 
228 
230 
233  int equilibrateMatrix();
234 
236 
239  int equilibrateRHS();
240 
242 
245  int applyRefinement();
246 
248 
251  int unequilibrateLHS();
252 
254 
262 
264 
265 
267  bool transpose() {return(transpose_);}
268 
270  bool factored() {return(factored_);}
271 
273  bool equilibratedA() {return(equilibratedA_);}
274 
276  bool equilibratedB() {return(equilibratedB_);}
277 
280 
283 
286 
288  bool solved() {return(solved_);}
289 
293 
295 
296 
299 
302 
305 
308 
310  OrdinalType numRows() const {return(M_);}
311 
313  OrdinalType numCols() const {return(N_);}
314 
316  OrdinalType numLower() const {return(KL_);}
317 
319  OrdinalType numUpper() const {return(KU_);}
320 
322  std::vector<OrdinalType> IPIV() const {return(IPIV_);}
323 
325  MagnitudeType ANORM() const {return(ANORM_);}
326 
328  MagnitudeType RCOND() const {return(RCOND_);}
329 
331 
333  MagnitudeType ROWCND() const {return(ROWCND_);}
334 
336 
338  MagnitudeType COLCND() const {return(COLCND_);}
339 
341  MagnitudeType AMAX() const {return(AMAX_);}
342 
344  std::vector<MagnitudeType> FERR() const {return(FERR_);}
345 
347  std::vector<MagnitudeType> BERR() const {return(BERR_);}
348 
350  std::vector<MagnitudeType> R() const {return(R_);}
351 
353  std::vector<MagnitudeType> C() const {return(C_);}
355 
357 
358  void Print(std::ostream& os) const;
361  protected:
362 
363  void allocateWORK() { LWORK_ = 3*N_; WORK_.resize( LWORK_ ); return;}
364  void resetMatrix();
365  void resetVectors();
366 
367 
373  bool factored_;
376  bool solved_;
380 
382 
383  OrdinalType M_;
384  OrdinalType N_;
385  OrdinalType KL_;
386  OrdinalType KU_;
387  OrdinalType Min_MN_;
388  OrdinalType LDA_;
389  OrdinalType LDAF_;
390  OrdinalType INFO_;
391  OrdinalType LWORK_;
392 
393  std::vector<OrdinalType> IPIV_;
394  std::vector<int> IWORK_;
395 
401 
406 
407  ScalarType * A_;
408  ScalarType * AF_;
409  std::vector<MagnitudeType> FERR_;
410  std::vector<MagnitudeType> BERR_;
411  std::vector<ScalarType> WORK_;
412  std::vector<MagnitudeType> R_;
413  std::vector<MagnitudeType> C_;
414 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
415  Eigen::PartialPivLU<EigenMatrix> lu_;
416 #endif
417 
418  private:
419 
420  // SerialBandDenseSolver copy constructor (put here because we don't want user access)
423 
424  };
425 
426  // Helper traits to distinguish work arrays for real and complex-valued datatypes.
427  using namespace details;
428 
429 //=============================================================================
430 
431 template<typename OrdinalType, typename ScalarType>
433  : CompObject(),
434  equilibrate_(false),
435  shouldEquilibrate_(false),
436  equilibratedA_(false),
437  equilibratedB_(false),
438  transpose_(false),
439  factored_(false),
440  estimateSolutionErrors_(false),
441  solutionErrorsEstimated_(false),
442  solved_(false),
443  reciprocalConditionEstimated_(false),
444  refineSolution_(false),
445  solutionRefined_(false),
446  TRANS_(Teuchos::NO_TRANS),
447  M_(0),
448  N_(0),
449  KL_(0),
450  KU_(0),
451  Min_MN_(0),
452  LDA_(0),
453  LDAF_(0),
454  INFO_(0),
455  LWORK_(0),
456  ANORM_(0.0),
457  RCOND_(ScalarTraits<MagnitudeType>::zero()),
458  ROWCND_(ScalarTraits<MagnitudeType>::zero()),
459  COLCND_(ScalarTraits<MagnitudeType>::zero()),
460  AMAX_(ScalarTraits<MagnitudeType>::zero()),
461  A_(NULL),
462  AF_(NULL)
463 {
464  resetMatrix();
465 }
466 //=============================================================================
467 
468 template<typename OrdinalType, typename ScalarType>
470 {}
471 
472 //=============================================================================
473 
474 template<typename OrdinalType, typename ScalarType>
476 {
477  LHS_ = Teuchos::null;
478  RHS_ = Teuchos::null;
479  reciprocalConditionEstimated_ = false;
480  solutionRefined_ = false;
481  solved_ = false;
482  solutionErrorsEstimated_ = false;
483  equilibratedB_ = false;
484 }
485 //=============================================================================
486 
487 template<typename OrdinalType, typename ScalarType>
489 {
490  resetVectors();
491  equilibratedA_ = false;
492  factored_ = false;
493  M_ = 0;
494  N_ = 0;
495  KL_ = 0;
496  KU_ = 0;
497  Min_MN_ = 0;
498  LDA_ = 0;
499  LDAF_ = 0;
504  A_ = 0;
505  AF_ = 0;
506  INFO_ = 0;
507  LWORK_ = 0;
508  R_.resize(0);
509  C_.resize(0);
510 }
511 //=============================================================================
512 
513 template<typename OrdinalType, typename ScalarType>
515 {
516 
517  OrdinalType m = AB->numRows();
518  OrdinalType n = AB->numCols();
519  OrdinalType kl = AB->lowerBandwidth();
520  OrdinalType ku = AB->upperBandwidth();
521 
522  // Check that the new matrix is consistent.
523  TEUCHOS_TEST_FOR_EXCEPTION(AB->values()==0, std::invalid_argument,
524  "SerialBandDenseSolver<T>::setMatrix: A is an empty SerialBandDenseMatrix<T>!");
525 
526  resetMatrix();
527  Matrix_ = AB;
528  Factor_ = AB;
529  M_ = m;
530  N_ = n;
531  KL_ = kl;
532  KU_ = ku-kl;
533  Min_MN_ = TEUCHOS_MIN(M_,N_);
534  LDA_ = AB->stride();
535  LDAF_ = LDA_;
536  A_ = AB->values();
537  AF_ = AB->values();
538 
539  return(0);
540 }
541 //=============================================================================
542 
543 template<typename OrdinalType, typename ScalarType>
546 {
547  // Check that these new vectors are consistent.
548  TEUCHOS_TEST_FOR_EXCEPTION(B->numRows()!=X->numRows() || B->numCols() != X->numCols(), std::invalid_argument,
549  "SerialBandDenseSolver<T>::setVectors: X and B are not the same size!");
550  TEUCHOS_TEST_FOR_EXCEPTION(B->values()==0, std::invalid_argument,
551  "SerialBandDenseSolver<T>::setVectors: B is an empty SerialDenseMatrix<T>!");
552  TEUCHOS_TEST_FOR_EXCEPTION(X->values()==0, std::invalid_argument,
553  "SerialBandDenseSolver<T>::setVectors: X is an empty SerialDenseMatrix<T>!");
554  TEUCHOS_TEST_FOR_EXCEPTION(B->stride()<1, std::invalid_argument,
555  "SerialBandDenseSolver<T>::setVectors: B has an invalid stride!");
556  TEUCHOS_TEST_FOR_EXCEPTION(X->stride()<1, std::invalid_argument,
557  "SerialBandDenseSolver<T>::setVectors: X has an invalid stride!");
558 
559  resetVectors();
560  LHS_ = X;
561  RHS_ = B;
562  return(0);
563 }
564 //=============================================================================
565 
566 template<typename OrdinalType, typename ScalarType>
568 {
569  estimateSolutionErrors_ = flag;
570 
571  // If the errors are estimated, this implies that the solution must be refined
572  refineSolution_ = refineSolution_ || flag;
573 }
574 //=============================================================================
575 
576 template<typename OrdinalType, typename ScalarType>
578 
579  if (factored()) return(0); // Already factored
580 
581  ANORM_ = Matrix_->normOne(); // Compute 1-Norm of A
582 
583  // If we want to refine the solution, then the factor must
584  // be stored separatedly from the original matrix
585 
586  if (A_ == AF_)
587  if (refineSolution_ ) {
588  Factor_ = rcp( new SerialBandDenseMatrix<OrdinalType,ScalarType>(*Matrix_) );
589  AF_ = Factor_->values();
590  LDAF_ = Factor_->stride();
591  }
592 
593  int ierr = 0;
594  if (equilibrate_) ierr = equilibrateMatrix();
595 
596  if (ierr!=0) return(ierr);
597 
598  if (IPIV_.size()==0) IPIV_.resize( N_ ); // Allocated Pivot vector if not already done.
599 
600  INFO_ = 0;
601 
602 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
603  EigenMatrixMap aMap( AF_+KL_, KL_+KU_+1, N_, EigenOuterStride(LDAF_) );
604  EigenMatrix fullMatrix(M_, N_);
605  for (OrdinalType j=0; j<N_; j++) {
606  for (OrdinalType i=TEUCHOS_MAX(0,j-KU_); i<=TEUCHOS_MIN(M_-1, j+KL_); i++) {
607  fullMatrix(i,j) = aMap(KU_+i-j, j);
608  }
609  }
610 
611  EigenPermutationMatrix p;
612  EigenOrdinalArray ind;
613  EigenOrdinalArrayMap ipivMap( &IPIV_[0], Min_MN_ );
614 
615  lu_.compute(fullMatrix);
616  p = lu_.permutationP();
617  ind = p.indices();
618 
619  for (OrdinalType i=0; i<ipivMap.innerSize(); i++) {
620  ipivMap(i) = ind(i);
621  }
622 
623 #else
624  this->GBTRF(M_, N_, KL_, KU_, AF_, LDAF_, &IPIV_[0], &INFO_);
625 #endif
626 
627  factored_ = true;
628 
629  return(INFO_);
630 }
631 
632 //=============================================================================
633 
634 template<typename OrdinalType, typename ScalarType>
636 
637  // If the user want the matrix to be equilibrated or wants a refined solution, we will
638  // call the X interface.
639  // Otherwise, if the matrix is already factored we will call the TRS interface.
640  // Otherwise, if the matrix is unfactored we will call the SV interface.
641 
642  int ierr = 0;
643  if (equilibrate_) {
644  ierr = equilibrateRHS();
645  }
646  if (ierr != 0) return(ierr); // Can't equilibrate B, so return.
647 
648  TEUCHOS_TEST_FOR_EXCEPTION( RHS_==Teuchos::null, std::invalid_argument,
649  "SerialBandDenseSolver<T>::solve: No right-hand side vector (RHS) has been set for the linear system!");
650  TEUCHOS_TEST_FOR_EXCEPTION( LHS_==Teuchos::null, std::invalid_argument,
651  "SerialBandDenseSolver<T>::solve: No solution vector (LHS) has been set for the linear system!");
652 
653  if (!factored()) factor(); // Matrix must be factored
654 
655  TEUCHOS_TEST_FOR_EXCEPTION( (equilibratedA_ && !equilibratedB_) || (!equilibratedA_ && equilibratedB_) ,
656  std::logic_error, "SerialBandDenseSolver<T>::solve: Matrix and vectors must be similarly scaled!");
657 
658  if (shouldEquilibrate() && !equilibratedA_)
659  std::cout << "WARNING! SerialBandDenseSolver<T>::solve: System should be equilibrated!" << std::endl;
660 
661  if (RHS_->values()!=LHS_->values()) {
662  (*LHS_) = (*RHS_); // Copy B to X if needed
663  }
664  INFO_ = 0;
665 
666 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
667  EigenMatrixMap rhsMap( RHS_->values(), RHS_->numRows(), RHS_->numCols(), EigenOuterStride(RHS_->stride()) );
668  EigenMatrixMap lhsMap( LHS_->values(), LHS_->numRows(), LHS_->numCols(), EigenOuterStride(LHS_->stride()) );
669  if ( TRANS_ == Teuchos::NO_TRANS ) {
670  lhsMap = lu_.solve(rhsMap);
671  } else if ( TRANS_ == Teuchos::TRANS ) {
672  lu_.matrixLU().template triangularView<Eigen::Upper>().transpose().solveInPlace(lhsMap);
673  lu_.matrixLU().template triangularView<Eigen::UnitLower>().transpose().solveInPlace(lhsMap);
674  EigenMatrix x = lu_.permutationP().transpose() * lhsMap;
675  lhsMap = x;
676  } else {
677  lu_.matrixLU().template triangularView<Eigen::Upper>().adjoint().solveInPlace(lhsMap);
678  lu_.matrixLU().template triangularView<Eigen::UnitLower>().adjoint().solveInPlace(lhsMap);
679  EigenMatrix x = lu_.permutationP().transpose() * lhsMap;
680  lhsMap = x;
681  }
682 #else
683  this->GBTRS(ETranspChar[TRANS_], N_, KL_, KU_, RHS_->numCols(), AF_, LDAF_, &IPIV_[0], LHS_->values(), LHS_->stride(), &INFO_);
684 #endif
685 
686  if (INFO_!=0) return(INFO_);
687  solved_ = true;
688 
689  int ierr1=0;
690  if (refineSolution_) ierr1 = applyRefinement();
691  if (ierr1!=0)
692  return(ierr1);
693 
694  if (equilibrate_) ierr1 = unequilibrateLHS();
695 
696  return(ierr1);
697 }
698 //=============================================================================
699 
700 template<typename OrdinalType, typename ScalarType>
702 {
703  TEUCHOS_TEST_FOR_EXCEPTION(!solved(), std::logic_error,
704  "SerialBandDenseSolver<T>::applyRefinement: Must have an existing solution!");
705  TEUCHOS_TEST_FOR_EXCEPTION(A_==AF_, std::logic_error,
706  "SerialBandDenseSolver<T>::applyRefinement: Cannot apply refinement if no original copy of A!");
707 
708 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
709  // Implement templated GERFS or use Eigen.
710  return (-1);
711 #else
712 
713  OrdinalType NRHS = RHS_->numCols();
714  FERR_.resize( NRHS );
715  BERR_.resize( NRHS );
716  allocateWORK();
717 
718  INFO_ = 0;
719  std::vector<typename details::lapack_traits<ScalarType>::iwork_type> GBRFS_WORK( N_ );
720  this->GBRFS(ETranspChar[TRANS_], N_, KL_, KU_, NRHS, A_+KL_, LDA_, AF_, LDAF_, &IPIV_[0],
721  RHS_->values(), RHS_->stride(), LHS_->values(), LHS_->stride(),
722  &FERR_[0], &BERR_[0], &WORK_[0], &GBRFS_WORK[0], &INFO_);
723 
724  solutionErrorsEstimated_ = true;
725  reciprocalConditionEstimated_ = true;
726  solutionRefined_ = true;
727 
728  return(INFO_);
729 #endif
730 }
731 
732 //=============================================================================
733 
734 template<typename OrdinalType, typename ScalarType>
736 {
737  if (R_.size()!=0) return(0); // Already computed
738 
739  R_.resize( M_ );
740  C_.resize( N_ );
741 
742  INFO_ = 0;
743  this->GBEQU (M_, N_, KL_, KU_, AF_+KL_, LDAF_, &R_[0], &C_[0], &ROWCND_, &COLCND_, &AMAX_, &INFO_);
744 
745  if (COLCND_<0.1*ScalarTraits<MagnitudeType>::one() ||
746  ROWCND_<0.1*ScalarTraits<MagnitudeType>::one() ||
748  shouldEquilibrate_ = true;
749 
750  return(INFO_);
751 }
752 
753 //=============================================================================
754 
755 template<typename OrdinalType, typename ScalarType>
757 {
758  OrdinalType i, j;
759  int ierr = 0;
760 
761  if (equilibratedA_) return(0); // Already done.
762  if (R_.size()==0) ierr = computeEquilibrateScaling(); // Compute R and C if needed.
763  if (ierr!=0) return(ierr); // If return value is not zero, then can't equilibrate.
764 
765  if (A_==AF_) {
766 
767  ScalarType * ptr;
768  for (j=0; j<N_; j++) {
769  ptr = A_ + KL_ + j*LDA_ + TEUCHOS_MAX(KU_-j, 0);
770  ScalarType s1 = C_[j];
771  for (i=TEUCHOS_MAX(0,j-KU_); i<=TEUCHOS_MIN(M_-1,j+KL_); i++) {
772  *ptr = *ptr*s1*R_[i];
773  ptr++;
774  }
775  }
776  } else {
777 
778  ScalarType * ptr;
779  ScalarType * ptrL;
780  ScalarType * ptrU;
781  for (j=0; j<N_; j++) {
782  ptr = A_ + KL_ + j*LDA_ + TEUCHOS_MAX(KU_-j, 0);
783  ScalarType s1 = C_[j];
784  for (i=TEUCHOS_MAX(0,j-KU_); i<=TEUCHOS_MIN(M_-1,j+KL_); i++) {
785  *ptr = *ptr*s1*R_[i];
786  ptr++;
787  }
788  }
789  for (j=0; j<N_; j++) {
790  ptrU = AF_ + j*LDAF_ + TEUCHOS_MAX(KL_+KU_-j, 0);
791  ptrL = AF_ + KL_ + KU_ + 1 + j*LDAF_;
792  ScalarType s1 = C_[j];
793  for (i=TEUCHOS_MAX(0,j-(KL_+KU_)); i<=TEUCHOS_MIN(M_-1,j); i++) {
794  *ptrU = *ptrU*s1*R_[i];
795  ptrU++;
796  }
797  for (i=TEUCHOS_MAX(0,j); i<=TEUCHOS_MIN(M_-1,j+KL_); i++) {
798  *ptrL = *ptrL*s1*R_[i];
799  ptrL++;
800  }
801  }
802  }
803 
804  equilibratedA_ = true;
805 
806  return(0);
807 }
808 
809 //=============================================================================
810 
811 template<typename OrdinalType, typename ScalarType>
813 {
814  OrdinalType i, j;
815  int ierr = 0;
816 
817  if (equilibratedB_) return(0); // Already done.
818  if (R_.size()==0) ierr = computeEquilibrateScaling(); // Compute R and C if needed.
819  if (ierr!=0) return(ierr); // Can't count on R and C being computed.
820 
821  MagnitudeType * R_tmp = &R_[0];
822  if (transpose_) R_tmp = &C_[0];
823 
824  OrdinalType LDB = RHS_->stride(), NRHS = RHS_->numCols();
825  ScalarType * B = RHS_->values();
826  ScalarType * ptr;
827  for (j=0; j<NRHS; j++) {
828  ptr = B + j*LDB;
829  for (i=0; i<M_; i++) {
830  *ptr = *ptr*R_tmp[i];
831  ptr++;
832  }
833  }
834 
835  equilibratedB_ = true;
836 
837  return(0);
838 }
839 
840 
841 //=============================================================================
842 
843 template<typename OrdinalType, typename ScalarType>
845 {
846  OrdinalType i, j;
847 
848  if (!equilibratedB_) return(0); // Nothing to do
849 
850  MagnitudeType * C_tmp = &C_[0];
851  if (transpose_) C_tmp = &R_[0];
852 
853  OrdinalType LDX = LHS_->stride(), NLHS = LHS_->numCols();
854  ScalarType * X = LHS_->values();
855  ScalarType * ptr;
856  for (j=0; j<NLHS; j++) {
857  ptr = X + j*LDX;
858  for (i=0; i<N_; i++) {
859  *ptr = *ptr*C_tmp[i];
860  ptr++;
861  }
862  }
863 
864  return(0);
865 }
866 
867 //=============================================================================
868 
869 template<typename OrdinalType, typename ScalarType>
871 {
872 #ifdef HAVE_TEUCHOSNUMERICS_EIGEN
873  // Implement templated GECON or use Eigen. Eigen currently doesn't have a condition estimation function.
874  return (-1);
875 #else
876 
877  if (reciprocalConditionEstimated()) {
878  Value = RCOND_;
879  return(0); // Already computed, just return it.
880  }
881 
882  if ( ANORM_<ScalarTraits<MagnitudeType>::zero() ) ANORM_ = Matrix_->normOne();
883 
884  int ierr = 0;
885  if (!factored()) ierr = factor(); // Need matrix factored.
886  if (ierr!=0) return(ierr);
887 
888  allocateWORK();
889 
890  // We will assume a one-norm condition number
891  INFO_ = 0;
892  std::vector<typename details::lapack_traits<ScalarType>::iwork_type> GBCON_WORK( N_ );
893  this->GBCON( '1', N_, KL_, KU_, AF_, LDAF_, &IPIV_[0], ANORM_, &RCOND_, &WORK_[0], &GBCON_WORK[0], &INFO_);
894  reciprocalConditionEstimated_ = true;
895  Value = RCOND_;
896 
897  return(INFO_);
898 #endif
899 }
900 //=============================================================================
901 
902 template<typename OrdinalType, typename ScalarType>
904 
905  if (Matrix_!=Teuchos::null) os << "Solver Matrix" << std::endl << *Matrix_ << std::endl;
906  if (Factor_!=Teuchos::null) os << "Solver Factored Matrix" << std::endl << *Factor_ << std::endl;
907  if (LHS_ !=Teuchos::null) os << "Solver LHS" << std::endl << *LHS_ << std::endl;
908  if (RHS_ !=Teuchos::null) os << "Solver RHS" << std::endl << *RHS_ << std::endl;
909 
910 }
911 
912 //=============================================================================
913 
914 
915 //=============================================================================
916 
917 
918 } // namespace Teuchos
919 
920 #endif /* _TEUCHOS_SERIALBANDDENSESOLVER_HPP_ */
bool equilibratedA()
Returns true if factor is equilibrated (factor available via AF() and LDAF()).
MagnitudeType RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed)...
RCP< SerialBandDenseMatrix< OrdinalType, ScalarType > > getFactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
void estimateSolutionErrors(bool flag)
Causes all solves to estimate the forward and backward solution error.
T magnitudeType
Mandatory typedef for result of magnitude.
Templated serial dense matrix class.
int applyRefinement()
Apply Iterative Refinement.
SerialBandDenseSolver & operator=(const SerialBandDenseSolver< OrdinalType, ScalarType > &Source)
void solveToRefinedSolution(bool flag)
Set whether or not to use iterative refinement to improve solutions to linear systems.
int setMatrix(const RCP< SerialBandDenseMatrix< OrdinalType, ScalarType > > &AB)
Sets the pointer for coefficient matrix.
Templated interface class to BLAS routines.
std::vector< MagnitudeType > R() const
Returns a pointer to the row scaling vector used for equilibration.
bool reciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
MagnitudeType ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed)...
int reciprocalConditionEstimate(MagnitudeType &Value)
Returns the reciprocal of the 1-norm condition number of the this matrix.
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > RHS_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Macro for throwing an exception with breakpointing to ease debugging.
Teuchos header file which uses auto-configuration information to include necessary C++ headers...
Templated BLAS wrapper.
std::vector< OrdinalType > IPIV() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
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 unequilibrateLHS()
Unscales the solution vectors if equilibration was used to solve the system.
int equilibrateMatrix()
Equilibrates the this matrix.
Object for storing data and providing functionality that is common to all computational classes...
This structure defines some basic traits for a scalar field type.
int computeEquilibrateScaling()
Computes the scaling vector S(i) = 1/sqrt(A(i,i)) of the this matrix.
bool solutionRefined()
Returns true if the current set of vectors has been refined.
Templated class for solving dense linear problems.
Functionality and data that is common to all computational classes.
Templated interface class to LAPACK routines.
Templated serial dense matrix class.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Deprecated.
bool solutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
OrdinalType numLower() const
Returns lower bandwidth of system.
The base Teuchos class.
std::vector< MagnitudeType > FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
This class creates and provides basic support for banded dense matrices of templated type...
bool equilibratedB()
Returns true if RHS is equilibrated (RHS available via B() and LDB()).
A class for representing and solving banded dense linear systems.
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > getLHS() const
Returns pointer to current LHS.
The Templated LAPACK Wrapper Class.
#define TEUCHOS_MAX(x, y)
void Print(std::ostream &os) const
Print service methods; defines behavior of ostream &lt;&lt; operator.
bool transpose()
Returns true if transpose of this matrix has and will be used.
OrdinalType numRows() const
Returns row dimension of system.
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
int solve()
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()..
RCP< SerialBandDenseMatrix< OrdinalType, ScalarType > > Factor_
MagnitudeType ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
RCP< SerialBandDenseMatrix< OrdinalType, ScalarType > > Matrix_
bool factored()
Returns true if matrix is factored (factor available via AF() and LDAF()).
bool solved()
Returns true if the current set of vectors has been solved.
int factor()
Computes the in-place LU factorization of the matrix using the LAPACK routine _GBTRF.
bool shouldEquilibrate()
Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system...
MagnitudeType COLCND() const
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
virtual ~SerialBandDenseSolver()
SerialBandDenseSolver destructor.
SerialBandDenseSolver()
Default constructor; matrix should be set using setMatrix(), LHS and RHS set with setVectors()...
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[]
void solveWithTransposeFlag(Teuchos::ETransp trans)
Smart reference counting pointer class for automatic garbage collection.
std::vector< MagnitudeType > BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
OrdinalType numUpper() const
Returns upper bandwidth of system.
std::vector< MagnitudeType > C() const
Returns a pointer to the column scale vector used for equilibration.
#define TEUCHOS_MIN(x, y)
RCP< SerialBandDenseMatrix< OrdinalType, ScalarType > > getMatrix() const
Returns pointer to current matrix.
int equilibrateRHS()
Equilibrates the current RHS.
OrdinalType numCols() const
Returns column dimension of system.
RCP< SerialDenseMatrix< OrdinalType, ScalarType > > LHS_
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 > > getRHS() const
Returns pointer to current RHS.
This class creates and provides basic support for dense rectangular matrix of templated type...