/usr/include/trilinos/Teuchos_SerialBandDenseSolver.hpp is in libtrilinos-teuchos-dev 12.4.2-2.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 | // @HEADER
// ***********************************************************************
//
// Teuchos: Common Tools Package
// Copyright (2004) Sandia Corporation
//
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
//
// ***********************************************************************
// @HEADER
#ifndef _TEUCHOS_SERIALBANDDENSESOLVER_HPP_
#define _TEUCHOS_SERIALBANDDENSESOLVER_HPP_
/// \file Teuchos_SerialBandDenseSolver.hpp
///
/// Declaration and definition of SerialBandDenseSolver,
/// a templated class for solving banded dense linear systems.
#include "Teuchos_CompObject.hpp"
#include "Teuchos_BLAS.hpp"
#include "Teuchos_LAPACK.hpp"
#include "Teuchos_RCP.hpp"
#include "Teuchos_ConfigDefs.hpp"
#include "Teuchos_SerialDenseMatrix.hpp"
#include "Teuchos_SerialBandDenseMatrix.hpp"
#include "Teuchos_SerialDenseSolver.hpp"
#include "Teuchos_ScalarTraits.hpp"
#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
#include "Eigen/Dense"
#endif
namespace Teuchos {
/*! \class SerialBandDenseSolver
\brief A class for representing and solving banded dense linear systems.
\tparam OrdinalType The index type used by the linear algebra implementation.
This should always be \c int.
\tparam ScalarType The type of entries in the matrix.
\section Teuchos_SerialBandDenseSolver_Intro Introduction
This class defines a banded dense matrix, which may have any number
of rows or columns (not necessarily equal). It's called "serial"
because the matrix lives in a single memory space. Thus, it's the
kind of matrix that one might give to the BLAS or LAPACK, not a
distributed matrix like one would give to ScaLAPACK.
This class also has methods for computing the (banded) LU
factorization of the matrix, and solving linear systems with the
matrix. We use instances of SerialDenseVector to represent the
right-hand side b or the solution x in the linear system \f$Ax=b\f$.
The instance of this class used store the banded matrix must contain
KL extra superdiagonals to store the L and U factors (see details
below).
Users have the option to do equilibration before factoring the matrix.
This may improve accuracy when solving ill-conditioned problems.
\section Teuchos_SerialBandDenseSolver_LAPACK SerialBandDenseSolver and LAPACK
Teuchos' LAPACK class wraps LAPACK's LU factorizations, including
the banded factorizations. It thus gives access to most of the same
functionality as SerialBandDenseSolver. The main difference is that
SerialBandDenseSolver offers a higher level of abstraction. It
hides the details of which LAPACK routines to call. Furthermore, if
you have built Teuchos with support for the third-party library
Eigen, SerialBandDenseSolver lets you solve linear systems for
<tt>ScalarType</tt> other than the four supported by the LAPACK
library.
\section Teuchos_SerialBandDenseSolver_Construct Constructing SerialBandDenseSolver objects
There is a single Teuchos::SerialBandDenseSolver constructor.
However, the matrix, right hand side and solution vectors must be
set prior to executing most methods in this class.
\subsection Teuchos_SerialBandDenseSolver_Construct_Setting Setting vectors used for linear solves
The matrix A, the left hand side X and the right hand side B (when
solving AX = B, for X), can be set by appropriate set methods. Each
of these three objects must be a SerialDenseMatrix or a
SerialDenseVector object. The set methods are as follows:
- setMatrix(): Sets the matrix
- setVectors() - Sets the left and right hand side vector(s)
\subsection Teuchos_SerialBandDenseSolver_Construct_Format Format of the matrix A
The SerialBandDenseMatrix must contain KL extra superdiagonals to store the L and U factors, where KL
is the upper bandwidth. Consider using the non-member conversion routines generalToBanded and bandedToGeneral if the
full SerialDenseMatrix is already in storage. However, it is more efficient simply to construct the
SerialBandDenseMatrix with the desired parameters and use the provided matrix access operators so
that the full rectangular matrix need not be stored. The conversion routine generalToBanded has a flag to store
the input Teuchos::SerialDenseMatrix in banded format with KL extra superdiagonals so this class can use it. Again,
it is more efficient to simply construct a Teuchos::SerialBandDenseMatrix object with KL extra superdiagonals than are
needed for the matrix data and fill the matrix using the matrix access operators.
See the documentation of Teuchos::SerialBandDenseMatrix for further details on the storage format.
\section Teuchos_SerialBandDenseSolver_Util Vector and Utility Functions
Once a Teuchos::SerialBandDenseSolver is constructed, several mathematical functions can be applied to
the object. Specifically:
<ul>
<li> Conversion between storage formats
<li> Factorizations
<li> Solves
<li> Condition estimates
<li> Equilibration
<li> Norms
</ul>
\section Teuchos_SerialBandDenseSolver_Strategies Strategies for Solving Linear Systems
In many cases, linear systems can be accurately solved by simply computing the LU factorization
of the matrix and then performing a forward back solve with a given set of right hand side vectors. However,
in some instances, the factorization may be very poorly conditioned and this simple approach may not work. In
these situations, equilibration and iterative refinement may improve the accuracy, or prevent a breakdown in
the factorization.
SerialBandDenseSolver will use equilibration with the factorization if, once the object
is constructed and \e before it is factored, you call the function factorWithEquilibration(true) to force
equilibration to be used. If you are uncertain if equilibration should be used, you may call the function
shouldEquilibrate() which will return true if equilibration could possibly help. shouldEquilibrate() uses
guidelines specified in the LAPACK User Guide, namely if SCOND < 0.1 and AMAX < Underflow or AMAX > Overflow, to
determine if equilibration \e might be useful.
SerialBandDenseSolver will use iterative refinement after a forward/back solve if you call
solveToRefinedSolution(true). It will also compute forward and backward error estimates if you call
estimateSolutionErrors(true). Access to the forward (back) error estimates is available via FERR() (BERR()).
Examples using SerialBandDenseSolver can be found in the Teuchos test directories.
*/
template<typename OrdinalType, typename ScalarType>
class SerialBandDenseSolver : public CompObject, public Object, public BLAS<OrdinalType, ScalarType>,
public LAPACK<OrdinalType, ScalarType>
{
public:
typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
typedef typename Eigen::Matrix<ScalarType,Eigen::Dynamic,Eigen::Dynamic,Eigen::ColMajor> EigenMatrix;
typedef typename Eigen::internal::BandMatrix<ScalarType,Eigen::Dynamic,Eigen::Dynamic,Eigen::ColMajor> EigenBandMatrix;
typedef typename Eigen::Matrix<ScalarType,Eigen::Dynamic,1> EigenVector;
typedef typename Eigen::InnerStride<Eigen::Dynamic> EigenInnerStride;
typedef typename Eigen::OuterStride<Eigen::Dynamic> EigenOuterStride;
typedef typename Eigen::Map<EigenVector,0,EigenInnerStride > EigenVectorMap;
typedef typename Eigen::Map<const EigenVector,0,EigenInnerStride > EigenConstVectorMap;
typedef typename Eigen::Map<EigenMatrix,0,EigenOuterStride > EigenMatrixMap;
typedef typename Eigen::Map<const EigenMatrix,0,EigenOuterStride > EigenConstMatrixMap;
typedef typename Eigen::PermutationMatrix<Eigen::Dynamic,Eigen::Dynamic,OrdinalType> EigenPermutationMatrix;
typedef typename Eigen::Array<OrdinalType,Eigen::Dynamic,1> EigenOrdinalArray;
typedef typename Eigen::Map<EigenOrdinalArray> EigenOrdinalArrayMap;
typedef typename Eigen::Array<ScalarType,Eigen::Dynamic,1> EigenScalarArray;
typedef typename Eigen::Map<EigenScalarArray> EigenScalarArrayMap;
#endif
//! @name Constructor/Destructor Methods
//@{
//! Default constructor; matrix should be set using setMatrix(), LHS and RHS set with setVectors().
SerialBandDenseSolver();
//! SerialBandDenseSolver destructor.
virtual ~SerialBandDenseSolver();
//@}
//! @name Set Methods
//@{
//! Sets the pointer for coefficient matrix
int setMatrix(const RCP<SerialBandDenseMatrix<OrdinalType, ScalarType> >& AB);
//! Sets the pointers for left and right hand side vector(s).
/*! Row dimension of X must match column dimension of matrix A, row dimension of B
must match row dimension of A. X and B must have the same dimensions.
*/
int setVectors(const RCP<SerialDenseMatrix<OrdinalType, ScalarType> >& X,
const RCP<SerialDenseMatrix<OrdinalType, ScalarType> >& B);
//@}
//! @name Strategy Modifying Methods
//@{
/// Set whether or not to equilibrate just before the matrix factorization.
/// This function must be called before the factorization is performed.
void factorWithEquilibration(bool flag) {equilibrate_ = flag; return;}
/// If \c flag is true, causes all subsequent function calls to work with the transpose of \e this matrix, otherwise not.
///
/// \note This interface will not work correctly for complex-valued linear systems, use solveWithTransposeFlag().
void solveWithTranspose(bool flag) {transpose_ = flag; if (flag) TRANS_ = Teuchos::TRANS; else TRANS_ = Teuchos::NO_TRANS; return;}
/// All subsequent function calls will work with the transpose-type set by this method (\c Teuchos::NO_TRANS, \c Teuchos::TRANS, and \c Teuchos::CONJ_TRANS).
/// \note This interface will allow correct behavior for complex-valued linear systems, solveWithTranspose() will not.
void solveWithTransposeFlag(Teuchos::ETransp trans) {TRANS_ = trans; if (trans != Teuchos::NO_TRANS) { transpose_ = true; } }
//! Set whether or not to use iterative refinement to improve solutions to linear systems.
void solveToRefinedSolution(bool flag) {refineSolution_ = flag; return;}
//! Causes all solves to estimate the forward and backward solution error.
/*! Error estimates will be in the arrays FERR and BERR, resp, after the solve step is complete.
These arrays are accessible via the FERR() and BERR() access functions.
*/
void estimateSolutionErrors(bool flag);
//@}
//! @name Factor/Solve Methods
//@{
//! Computes the in-place LU factorization of the matrix using the LAPACK routine \e _GBTRF.
/*!
\return Integer error code, set to 0 if successful.
*/
int factor();
//! Computes the solution X to AX = B for the \e this matrix and the B provided to SetVectors()..
/*!
\return Integer error code, set to 0 if successful.
*/
int solve();
//! Computes the scaling vector S(i) = 1/sqrt(A(i,i)) of the \e this matrix.
/*!
\return Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.
*/
int computeEquilibrateScaling();
//! Equilibrates the \e this matrix.
/*!
\return Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.
*/
int equilibrateMatrix();
//! Equilibrates the current RHS.
/*!
\return Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.
*/
int equilibrateRHS();
//! Apply Iterative Refinement.
/*!
\return Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.
*/
int applyRefinement();
//! Unscales the solution vectors if equilibration was used to solve the system.
/*!
\return Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.
*/
int unequilibrateLHS();
//! Returns the reciprocal of the 1-norm condition number of the \e this matrix.
/*!
\param Value Out
On return contains the reciprocal of the 1-norm condition number of the \e this matrix.
\return Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.
*/
int reciprocalConditionEstimate(MagnitudeType & Value);
//@}
//! @name Query methods
//@{
//! Returns true if transpose of \e this matrix has and will be used.
bool transpose() {return(transpose_);}
//! Returns true if matrix is factored (factor available via AF() and LDAF()).
bool factored() {return(factored_);}
//! Returns true if factor is equilibrated (factor available via AF() and LDAF()).
bool equilibratedA() {return(equilibratedA_);}
//! Returns true if RHS is equilibrated (RHS available via B() and LDB()).
bool equilibratedB() {return(equilibratedB_);}
//! Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system.
bool shouldEquilibrate() {computeEquilibrateScaling(); return(shouldEquilibrate_);}
//! Returns true if forward and backward error estimated have been computed (available via FERR() and BERR()).
bool solutionErrorsEstimated() {return(solutionErrorsEstimated_);}
//! Returns true if the condition number of the \e this matrix has been computed (value available via ReciprocalConditionEstimate()).
bool reciprocalConditionEstimated() {return(reciprocalConditionEstimated_);}
//! Returns true if the current set of vectors has been solved.
bool solved() {return(solved_);}
//! Returns true if the current set of vectors has been refined.
bool solutionRefined() {return(solutionRefined_);}
//@}
//! @name Data Accessor methods
//@{
//! Returns pointer to current matrix.
RCP<SerialBandDenseMatrix<OrdinalType, ScalarType> > getMatrix() const {return(Matrix_);}
//! Returns pointer to factored matrix (assuming factorization has been performed).
RCP<SerialBandDenseMatrix<OrdinalType, ScalarType> > getFactoredMatrix() const {return(Factor_);}
//! Returns pointer to current LHS.
RCP<SerialDenseMatrix<OrdinalType, ScalarType> > getLHS() const {return(LHS_);}
//! Returns pointer to current RHS.
RCP<SerialDenseMatrix<OrdinalType, ScalarType> > getRHS() const {return(RHS_);}
//! Returns row dimension of system.
OrdinalType numRows() const {return(M_);}
//! Returns column dimension of system.
OrdinalType numCols() const {return(N_);}
//! Returns lower bandwidth of system.
OrdinalType numLower() const {return(KL_);}
//! Returns upper bandwidth of system.
OrdinalType numUpper() const {return(KU_);}
//! Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
std::vector<OrdinalType> IPIV() const {return(IPIV_);}
//! Returns the 1-Norm of the \e this matrix (returns -1 if not yet computed).
MagnitudeType ANORM() const {return(ANORM_);}
//! Returns the reciprocal of the condition number of the \e this matrix (returns -1 if not yet computed).
MagnitudeType RCOND() const {return(RCOND_);}
//! Ratio of smallest to largest row scale factors for the \e this matrix (returns -1 if not yet computed).
/*! If ROWCND() is >= 0.1 and AMAX() is not close to overflow or underflow, then equilibration is not needed.
*/
MagnitudeType ROWCND() const {return(ROWCND_);}
//! Ratio of smallest to largest column scale factors for the \e this matrix (returns -1 if not yet computed).
/*! If COLCND() is >= 0.1 then equilibration is not needed.
*/
MagnitudeType COLCND() const {return(COLCND_);}
//! Returns the absolute value of the largest entry of the \e this matrix (returns -1 if not yet computed).
MagnitudeType AMAX() const {return(AMAX_);}
//! Returns a pointer to the forward error estimates computed by LAPACK.
std::vector<MagnitudeType> FERR() const {return(FERR_);}
//! Returns a pointer to the backward error estimates computed by LAPACK.
std::vector<MagnitudeType> BERR() const {return(BERR_);}
//! Returns a pointer to the row scaling vector used for equilibration.
std::vector<MagnitudeType> R() const {return(R_);}
//! Returns a pointer to the column scale vector used for equilibration.
std::vector<MagnitudeType> C() const {return(C_);}
//@}
//! @name I/O methods
//@{
//! Print service methods; defines behavior of ostream << operator.
void Print(std::ostream& os) const;
//@}
protected:
void allocateWORK() { LWORK_ = 3*N_; WORK_.resize( LWORK_ ); return;}
void resetMatrix();
void resetVectors();
bool equilibrate_;
bool shouldEquilibrate_;
bool equilibratedA_;
bool equilibratedB_;
bool transpose_;
bool factored_;
bool estimateSolutionErrors_;
bool solutionErrorsEstimated_;
bool solved_;
bool reciprocalConditionEstimated_;
bool refineSolution_;
bool solutionRefined_;
Teuchos::ETransp TRANS_;
OrdinalType M_;
OrdinalType N_;
OrdinalType KL_;
OrdinalType KU_;
OrdinalType Min_MN_;
OrdinalType LDA_;
OrdinalType LDAF_;
OrdinalType INFO_;
OrdinalType LWORK_;
std::vector<OrdinalType> IPIV_;
std::vector<int> IWORK_;
MagnitudeType ANORM_;
MagnitudeType RCOND_;
MagnitudeType ROWCND_;
MagnitudeType COLCND_;
MagnitudeType AMAX_;
RCP<SerialBandDenseMatrix<OrdinalType, ScalarType> > Matrix_;
RCP<SerialDenseMatrix<OrdinalType, ScalarType> > LHS_;
RCP<SerialDenseMatrix<OrdinalType, ScalarType> > RHS_;
RCP<SerialBandDenseMatrix<OrdinalType, ScalarType> > Factor_;
ScalarType * A_;
ScalarType * AF_;
std::vector<MagnitudeType> FERR_;
std::vector<MagnitudeType> BERR_;
std::vector<ScalarType> WORK_;
std::vector<MagnitudeType> R_;
std::vector<MagnitudeType> C_;
#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
Eigen::PartialPivLU<EigenMatrix> lu_;
#endif
private:
// SerialBandDenseSolver copy constructor (put here because we don't want user access)
SerialBandDenseSolver(const SerialBandDenseSolver<OrdinalType, ScalarType>& Source);
SerialBandDenseSolver & operator=(const SerialBandDenseSolver<OrdinalType, ScalarType>& Source);
};
// Helper traits to distinguish work arrays for real and complex-valued datatypes.
using namespace details;
//=============================================================================
template<typename OrdinalType, typename ScalarType>
SerialBandDenseSolver<OrdinalType,ScalarType>::SerialBandDenseSolver()
: CompObject(),
equilibrate_(false),
shouldEquilibrate_(false),
equilibratedA_(false),
equilibratedB_(false),
transpose_(false),
factored_(false),
estimateSolutionErrors_(false),
solutionErrorsEstimated_(false),
solved_(false),
reciprocalConditionEstimated_(false),
refineSolution_(false),
solutionRefined_(false),
TRANS_(Teuchos::NO_TRANS),
M_(0),
N_(0),
KL_(0),
KU_(0),
Min_MN_(0),
LDA_(0),
LDAF_(0),
INFO_(0),
LWORK_(0),
ANORM_(0.0),
RCOND_(ScalarTraits<MagnitudeType>::zero()),
ROWCND_(ScalarTraits<MagnitudeType>::zero()),
COLCND_(ScalarTraits<MagnitudeType>::zero()),
AMAX_(ScalarTraits<MagnitudeType>::zero()),
A_(NULL),
AF_(NULL)
{
resetMatrix();
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
SerialBandDenseSolver<OrdinalType,ScalarType>::~SerialBandDenseSolver()
{}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
void SerialBandDenseSolver<OrdinalType,ScalarType>::resetVectors()
{
LHS_ = Teuchos::null;
RHS_ = Teuchos::null;
reciprocalConditionEstimated_ = false;
solutionRefined_ = false;
solved_ = false;
solutionErrorsEstimated_ = false;
equilibratedB_ = false;
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
void SerialBandDenseSolver<OrdinalType,ScalarType>::resetMatrix()
{
resetVectors();
equilibratedA_ = false;
factored_ = false;
M_ = 0;
N_ = 0;
KL_ = 0;
KU_ = 0;
Min_MN_ = 0;
LDA_ = 0;
LDAF_ = 0;
RCOND_ = -ScalarTraits<MagnitudeType>::one();
ROWCND_ = -ScalarTraits<MagnitudeType>::one();
COLCND_ = -ScalarTraits<MagnitudeType>::one();
AMAX_ = -ScalarTraits<MagnitudeType>::one();
A_ = 0;
AF_ = 0;
INFO_ = 0;
LWORK_ = 0;
R_.resize(0);
C_.resize(0);
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
int SerialBandDenseSolver<OrdinalType,ScalarType>::setMatrix(const RCP<SerialBandDenseMatrix<OrdinalType,ScalarType> >& AB)
{
OrdinalType m = AB->numRows();
OrdinalType n = AB->numCols();
OrdinalType kl = AB->lowerBandwidth();
OrdinalType ku = AB->upperBandwidth();
// Check that the new matrix is consistent.
TEUCHOS_TEST_FOR_EXCEPTION(AB->values()==0, std::invalid_argument,
"SerialBandDenseSolver<T>::setMatrix: A is an empty SerialBandDenseMatrix<T>!");
resetMatrix();
Matrix_ = AB;
Factor_ = AB;
M_ = m;
N_ = n;
KL_ = kl;
KU_ = ku-kl;
Min_MN_ = TEUCHOS_MIN(M_,N_);
LDA_ = AB->stride();
LDAF_ = LDA_;
A_ = AB->values();
AF_ = AB->values();
return(0);
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
int SerialBandDenseSolver<OrdinalType,ScalarType>::setVectors(const RCP<SerialDenseMatrix<OrdinalType,ScalarType> >& X,
const RCP<SerialDenseMatrix<OrdinalType,ScalarType> >& B)
{
// Check that these new vectors are consistent.
TEUCHOS_TEST_FOR_EXCEPTION(B->numRows()!=X->numRows() || B->numCols() != X->numCols(), std::invalid_argument,
"SerialBandDenseSolver<T>::setVectors: X and B are not the same size!");
TEUCHOS_TEST_FOR_EXCEPTION(B->values()==0, std::invalid_argument,
"SerialBandDenseSolver<T>::setVectors: B is an empty SerialDenseMatrix<T>!");
TEUCHOS_TEST_FOR_EXCEPTION(X->values()==0, std::invalid_argument,
"SerialBandDenseSolver<T>::setVectors: X is an empty SerialDenseMatrix<T>!");
TEUCHOS_TEST_FOR_EXCEPTION(B->stride()<1, std::invalid_argument,
"SerialBandDenseSolver<T>::setVectors: B has an invalid stride!");
TEUCHOS_TEST_FOR_EXCEPTION(X->stride()<1, std::invalid_argument,
"SerialBandDenseSolver<T>::setVectors: X has an invalid stride!");
resetVectors();
LHS_ = X;
RHS_ = B;
return(0);
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
void SerialBandDenseSolver<OrdinalType,ScalarType>::estimateSolutionErrors(bool flag)
{
estimateSolutionErrors_ = flag;
// If the errors are estimated, this implies that the solution must be refined
refineSolution_ = refineSolution_ || flag;
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
int SerialBandDenseSolver<OrdinalType,ScalarType>::factor() {
if (factored()) return(0); // Already factored
ANORM_ = Matrix_->normOne(); // Compute 1-Norm of A
// If we want to refine the solution, then the factor must
// be stored separatedly from the original matrix
if (A_ == AF_)
if (refineSolution_ ) {
Factor_ = rcp( new SerialBandDenseMatrix<OrdinalType,ScalarType>(*Matrix_) );
AF_ = Factor_->values();
LDAF_ = Factor_->stride();
}
int ierr = 0;
if (equilibrate_) ierr = equilibrateMatrix();
if (ierr!=0) return(ierr);
if (IPIV_.size()==0) IPIV_.resize( N_ ); // Allocated Pivot vector if not already done.
INFO_ = 0;
#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
EigenMatrixMap aMap( AF_+KL_, KL_+KU_+1, N_, EigenOuterStride(LDAF_) );
EigenMatrix fullMatrix(M_, N_);
for (OrdinalType j=0; j<N_; j++) {
for (OrdinalType i=TEUCHOS_MAX(0,j-KU_); i<=TEUCHOS_MIN(M_-1, j+KL_); i++) {
fullMatrix(i,j) = aMap(KU_+i-j, j);
}
}
EigenPermutationMatrix p;
EigenOrdinalArray ind;
EigenOrdinalArrayMap ipivMap( &IPIV_[0], Min_MN_ );
lu_.compute(fullMatrix);
p = lu_.permutationP();
ind = p.indices();
for (OrdinalType i=0; i<ipivMap.innerSize(); i++) {
ipivMap(i) = ind(i);
}
#else
this->GBTRF(M_, N_, KL_, KU_, AF_, LDAF_, &IPIV_[0], &INFO_);
#endif
factored_ = true;
return(INFO_);
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
int SerialBandDenseSolver<OrdinalType,ScalarType>::solve() {
// If the user want the matrix to be equilibrated or wants a refined solution, we will
// call the X interface.
// Otherwise, if the matrix is already factored we will call the TRS interface.
// Otherwise, if the matrix is unfactored we will call the SV interface.
int ierr = 0;
if (equilibrate_) {
ierr = equilibrateRHS();
equilibratedB_ = true;
}
if (ierr != 0) return(ierr); // Can't equilibrate B, so return.
TEUCHOS_TEST_FOR_EXCEPTION( (equilibratedA_ && !equilibratedB_) || (!equilibratedA_ && equilibratedB_) ,
std::logic_error, "SerialBandDenseSolver<T>::solve: Matrix and vectors must be similarly scaled!");
TEUCHOS_TEST_FOR_EXCEPTION( RHS_==Teuchos::null, std::invalid_argument,
"SerialBandDenseSolver<T>::solve: No right-hand side vector (RHS) has been set for the linear system!");
TEUCHOS_TEST_FOR_EXCEPTION( LHS_==Teuchos::null, std::invalid_argument,
"SerialBandDenseSolver<T>::solve: No solution vector (LHS) has been set for the linear system!");
if (shouldEquilibrate() && !equilibratedA_)
std::cout << "WARNING! SerialBandDenseSolver<T>::solve: System should be equilibrated!" << std::endl;
if (!factored()) factor(); // Matrix must be factored
if (RHS_->values()!=LHS_->values()) {
(*LHS_) = (*RHS_); // Copy B to X if needed
}
INFO_ = 0;
#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
EigenMatrixMap rhsMap( RHS_->values(), RHS_->numRows(), RHS_->numCols(), EigenOuterStride(RHS_->stride()) );
EigenMatrixMap lhsMap( LHS_->values(), LHS_->numRows(), LHS_->numCols(), EigenOuterStride(LHS_->stride()) );
if ( TRANS_ == Teuchos::NO_TRANS ) {
lhsMap = lu_.solve(rhsMap);
} else if ( TRANS_ == Teuchos::TRANS ) {
lu_.matrixLU().template triangularView<Eigen::Upper>().transpose().solveInPlace(lhsMap);
lu_.matrixLU().template triangularView<Eigen::UnitLower>().transpose().solveInPlace(lhsMap);
EigenMatrix x = lu_.permutationP().transpose() * lhsMap;
lhsMap = x;
} else {
lu_.matrixLU().template triangularView<Eigen::Upper>().adjoint().solveInPlace(lhsMap);
lu_.matrixLU().template triangularView<Eigen::UnitLower>().adjoint().solveInPlace(lhsMap);
EigenMatrix x = lu_.permutationP().transpose() * lhsMap;
lhsMap = x;
}
#else
this->GBTRS(ETranspChar[TRANS_], N_, KL_, KU_, RHS_->numCols(), AF_, LDAF_, &IPIV_[0], LHS_->values(), LHS_->stride(), &INFO_);
#endif
if (INFO_!=0) return(INFO_);
solved_ = true;
int ierr1=0;
if (refineSolution_) ierr1 = applyRefinement();
if (ierr1!=0)
return(ierr1);
if (equilibrate_) ierr1 = unequilibrateLHS();
return(ierr1);
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
int SerialBandDenseSolver<OrdinalType,ScalarType>::applyRefinement()
{
TEUCHOS_TEST_FOR_EXCEPTION(!solved(), std::logic_error,
"SerialBandDenseSolver<T>::applyRefinement: Must have an existing solution!");
TEUCHOS_TEST_FOR_EXCEPTION(A_==AF_, std::logic_error,
"SerialBandDenseSolver<T>::applyRefinement: Cannot apply refinement if no original copy of A!");
#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
// Implement templated GERFS or use Eigen.
return (-1);
#else
OrdinalType NRHS = RHS_->numCols();
FERR_.resize( NRHS );
BERR_.resize( NRHS );
allocateWORK();
INFO_ = 0;
std::vector<typename details::lapack_traits<ScalarType>::iwork_type> GBRFS_WORK( N_ );
this->GBRFS(ETranspChar[TRANS_], N_, KL_, KU_, NRHS, A_+KL_, LDA_, AF_, LDAF_, &IPIV_[0],
RHS_->values(), RHS_->stride(), LHS_->values(), LHS_->stride(),
&FERR_[0], &BERR_[0], &WORK_[0], &GBRFS_WORK[0], &INFO_);
solutionErrorsEstimated_ = true;
reciprocalConditionEstimated_ = true;
solutionRefined_ = true;
return(INFO_);
#endif
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
int SerialBandDenseSolver<OrdinalType,ScalarType>::computeEquilibrateScaling()
{
if (R_.size()!=0) return(0); // Already computed
R_.resize( M_ );
C_.resize( N_ );
INFO_ = 0;
this->GBEQU (M_, N_, KL_, KU_, AF_+KL_, LDAF_, &R_[0], &C_[0], &ROWCND_, &COLCND_, &AMAX_, &INFO_);
if (COLCND_<0.1*ScalarTraits<MagnitudeType>::one() ||
ROWCND_<0.1*ScalarTraits<MagnitudeType>::one() ||
AMAX_ < ScalarTraits<ScalarType>::rmin() || AMAX_ > ScalarTraits<ScalarType>::rmax() )
shouldEquilibrate_ = true;
return(INFO_);
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
int SerialBandDenseSolver<OrdinalType,ScalarType>::equilibrateMatrix()
{
OrdinalType i, j;
int ierr = 0;
if (equilibratedA_) return(0); // Already done.
if (R_.size()==0) ierr = computeEquilibrateScaling(); // Compute R and C if needed.
if (ierr!=0) return(ierr); // If return value is not zero, then can't equilibrate.
if (A_==AF_) {
ScalarType * ptr;
for (j=0; j<N_; j++) {
ptr = A_ + KL_ + j*LDA_ + TEUCHOS_MAX(KU_-j, 0);
ScalarType s1 = C_[j];
for (i=TEUCHOS_MAX(0,j-KU_); i<=TEUCHOS_MIN(M_-1,j+KL_); i++) {
*ptr = *ptr*s1*R_[i];
ptr++;
}
}
} else {
ScalarType * ptr;
ScalarType * ptrL;
ScalarType * ptrU;
for (j=0; j<N_; j++) {
ptr = A_ + KL_ + j*LDA_ + TEUCHOS_MAX(KU_-j, 0);
ScalarType s1 = C_[j];
for (i=TEUCHOS_MAX(0,j-KU_); i<=TEUCHOS_MIN(M_-1,j+KL_); i++) {
*ptr = *ptr*s1*R_[i];
ptr++;
}
}
for (j=0; j<N_; j++) {
ptrU = AF_ + j*LDAF_ + TEUCHOS_MAX(KL_+KU_-j, 0);
ptrL = AF_ + KL_ + KU_ + 1 + j*LDAF_;
ScalarType s1 = C_[j];
for (i=TEUCHOS_MAX(0,j-(KL_+KU_)); i<=TEUCHOS_MIN(M_-1,j); i++) {
*ptrU = *ptrU*s1*R_[i];
ptrU++;
}
for (i=TEUCHOS_MAX(0,j); i<=TEUCHOS_MIN(M_-1,j+KL_); i++) {
*ptrL = *ptrL*s1*R_[i];
ptrL++;
}
}
}
equilibratedA_ = true;
return(0);
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
int SerialBandDenseSolver<OrdinalType,ScalarType>::equilibrateRHS()
{
OrdinalType i, j;
int ierr = 0;
if (equilibratedB_) return(0); // Already done.
if (R_.size()==0) ierr = computeEquilibrateScaling(); // Compute R and C if needed.
if (ierr!=0) return(ierr); // Can't count on R and C being computed.
MagnitudeType * R_tmp = &R_[0];
if (transpose_) R_tmp = &C_[0];
OrdinalType LDB = RHS_->stride(), NRHS = RHS_->numCols();
ScalarType * B = RHS_->values();
ScalarType * ptr;
for (j=0; j<NRHS; j++) {
ptr = B + j*LDB;
for (i=0; i<M_; i++) {
*ptr = *ptr*R_tmp[i];
ptr++;
}
}
equilibratedB_ = true;
return(0);
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
int SerialBandDenseSolver<OrdinalType,ScalarType>::unequilibrateLHS()
{
OrdinalType i, j;
if (!equilibratedB_) return(0); // Nothing to do
MagnitudeType * C_tmp = &C_[0];
if (transpose_) C_tmp = &R_[0];
OrdinalType LDX = LHS_->stride(), NLHS = LHS_->numCols();
ScalarType * X = LHS_->values();
ScalarType * ptr;
for (j=0; j<NLHS; j++) {
ptr = X + j*LDX;
for (i=0; i<N_; i++) {
*ptr = *ptr*C_tmp[i];
ptr++;
}
}
return(0);
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
int SerialBandDenseSolver<OrdinalType,ScalarType>::reciprocalConditionEstimate(MagnitudeType & Value)
{
#ifdef HAVE_TEUCHOSNUMERICS_EIGEN
// Implement templated GECON or use Eigen. Eigen currently doesn't have a condition estimation function.
return (-1);
#else
if (reciprocalConditionEstimated()) {
Value = RCOND_;
return(0); // Already computed, just return it.
}
if ( ANORM_<ScalarTraits<MagnitudeType>::zero() ) ANORM_ = Matrix_->normOne();
int ierr = 0;
if (!factored()) ierr = factor(); // Need matrix factored.
if (ierr!=0) return(ierr);
allocateWORK();
// We will assume a one-norm condition number
INFO_ = 0;
std::vector<typename details::lapack_traits<ScalarType>::iwork_type> GBCON_WORK( N_ );
this->GBCON( '1', N_, KL_, KU_, AF_, LDAF_, &IPIV_[0], ANORM_, &RCOND_, &WORK_[0], &GBCON_WORK[0], &INFO_);
reciprocalConditionEstimated_ = true;
Value = RCOND_;
return(INFO_);
#endif
}
//=============================================================================
template<typename OrdinalType, typename ScalarType>
void SerialBandDenseSolver<OrdinalType,ScalarType>::Print(std::ostream& os) const {
if (Matrix_!=Teuchos::null) os << "Solver Matrix" << std::endl << *Matrix_ << std::endl;
if (Factor_!=Teuchos::null) os << "Solver Factored Matrix" << std::endl << *Factor_ << std::endl;
if (LHS_ !=Teuchos::null) os << "Solver LHS" << std::endl << *LHS_ << std::endl;
if (RHS_ !=Teuchos::null) os << "Solver RHS" << std::endl << *RHS_ << std::endl;
}
//=============================================================================
//=============================================================================
} // namespace Teuchos
#endif /* _TEUCHOS_SERIALBANDDENSESOLVER_HPP_ */
|