Logo Search packages:      
Sourcecode: eigen3 version File versions  Download package

BiCGSTAB.h
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2011 Gael Guennebaud <gael.guennebaud@inria.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.

#ifndef EIGEN_BICGSTAB_H
#define EIGEN_BICGSTAB_H

namespace internal {

/** \internal Low-level bi conjugate gradient stabilized algorithm
  * \param mat The matrix A
  * \param rhs The right hand side vector b
  * \param x On input and initial solution, on output the computed solution.
  * \param precond A preconditioner being able to efficiently solve for an
  *                approximation of Ax=b (regardless of b)
  * \param iters On input the max number of iteration, on output the number of performed iterations.
  * \param tol_error On input the tolerance error, on output an estimation of the relative error.
  */
template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
00040 void bicgstab(const MatrixType& mat, const Rhs& rhs, Dest& x,
              const Preconditioner& precond, int& iters,
              typename Dest::RealScalar& tol_error)
{
  using std::sqrt;
  using std::abs;
  typedef typename Dest::RealScalar RealScalar;
  typedef typename Dest::Scalar Scalar;
  typedef Matrix<Scalar,Dynamic,1> VectorType;
  
  RealScalar tol = tol_error;
  int maxIters = iters;

  int n = mat.cols();
  VectorType r  = rhs - mat * x;
  VectorType r0 = r;
  RealScalar r0_sqnorm = r0.squaredNorm();
  Scalar rho    = 1;
  Scalar alpha  = 1;
  Scalar w      = 1;
  
  VectorType v = VectorType::Zero(n), p = VectorType::Zero(n);
  VectorType y(n),  z(n);
  VectorType kt(n), ks(n);

  VectorType s(n), t(n);

  RealScalar tol2 = tol*tol;
  int i = 0;

  do
  {
    Scalar rho_old = rho;

    rho = r0.dot(r);
    Scalar beta = (rho/rho_old) * (alpha / w);
    p = r + beta * (p - w * v);
    
    y = precond.solve(p);
    v.noalias() = mat * y;

    alpha = rho / r0.dot(v);
    s = r - alpha * v;

    z = precond.solve(s);
    t.noalias() = mat * z;

    kt = precond.solve(t);
    ks = precond.solve(s);

    w = kt.dot(ks) / kt.squaredNorm();
    x += alpha * y + w * z;
    r = s - w * t;
    ++i;
  } while ( r.squaredNorm()/r0_sqnorm > tol2 && i<maxIters );
  
  tol_error = sqrt(r.squaredNorm()/r0_sqnorm);
  //tol_error = sqrt(abs(absNew / absInit));
  iters = i;
}

}

template< typename _MatrixType,
          typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
class BiCGSTAB;

namespace internal {

template< typename _MatrixType, typename _Preconditioner>
00110 struct traits<BiCGSTAB<_MatrixType,_Preconditioner> >
{
  typedef _MatrixType MatrixType;
  typedef _Preconditioner Preconditioner;
};

}

/** \ingroup IterativeLinearSolvers_Module
  * \brief A bi conjugate gradient stabilized solver for sparse square problems
  *
  * This class allows to solve for A.x = b sparse linear problems using a bi conjugate gradient
  * stabilized algorithm. The vectors x and b can be either dense or sparse.
  *
  * \tparam _MatrixType the type of the sparse matrix A, can be a dense or a sparse matrix.
  * \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
  *
  * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
  * and setTolerance() methods. The default are 1000 max iterations and NumTraits<Scalar>::epsilon()
  * for the tolerance.
  * 
  * This class can be used as the direct solver classes. Here is a typical usage example:
  * \code
  * int n = 10000;
  * VectorXd x(n), b(n);
  * SparseMatrix<double> A(n,n);
  * // fill A and b
  * BiCGSTAB<SparseMatrix<double> > solver;
  * solver(A);
  * x = solver.solve(b);
  * std::cout << "#iterations:     " << solver.iterations() << std::endl;
  * std::cout << "estimated error: " << solver.error()      << std::endl;
  * // update b, and solve again
  * x = solver.solve(b);
  * \endcode
  * 
  * By default the iterations start with x=0 as an initial guess of the solution.
  * One can control the start using the solveWithGuess() method. Here is a step by
  * step execution example starting with a random guess and printing the evolution
  * of the estimated error:
  * * \code
  * x = VectorXd::Random(n);
  * solver.setMaxIterations(1);
  * int i = 0;
  * do {
  *   x = solver.solveWithGuess(b,x);
  *   std::cout << i << " : " << solver.error() << std::endl;
  *   ++i;
  * } while (solver.info()!=Success && i<100);
  * \endcode
  * Note that such a step by step excution is slightly slower.
  * 
  * \sa class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
  */
template< typename _MatrixType, typename _Preconditioner>
00165 class BiCGSTAB : public IterativeSolverBase<BiCGSTAB<_MatrixType,_Preconditioner> >
{
  typedef IterativeSolverBase<BiCGSTAB> Base;
  using Base::mp_matrix;
  using Base::m_error;
  using Base::m_iterations;
  using Base::m_info;
  using Base::m_isInitialized;
public:
  typedef _MatrixType MatrixType;
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::Index Index;
  typedef typename MatrixType::RealScalar RealScalar;
  typedef _Preconditioner Preconditioner;

public:

  /** Default constructor. */
00183   BiCGSTAB() : Base() {}

  /** Initialize the solver with matrix \a A for further \c Ax=b solving.
    * 
    * This constructor is a shortcut for the default constructor followed
    * by a call to compute().
    * 
    * \warning this class stores a reference to the matrix A as well as some
    * precomputed values that depend on it. Therefore, if \a A is changed
    * this class becomes invalid. Call compute() to update it with the new
    * matrix A, or modify a copy of A.
    */
00195   BiCGSTAB(const MatrixType& A) : Base(A) {}

  ~BiCGSTAB() {}
  
  /** \returns the solution x of \f$ A x = b \f$ using the current decomposition of A
    * \a x0 as an initial solution.
    *
    * \sa compute()
    */
  template<typename Rhs,typename Guess>
  inline const internal::solve_retval_with_guess<BiCGSTAB, Rhs, Guess>
00206   solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
  {
    eigen_assert(m_isInitialized && "BiCGSTAB is not initialized.");
    eigen_assert(Base::rows()==b.rows()
              && "BiCGSTAB::solve(): invalid number of rows of the right hand side matrix b");
    return internal::solve_retval_with_guess
            <BiCGSTAB, Rhs, Guess>(*this, b.derived(), x0);
  }
  
  /** \internal */
  template<typename Rhs,typename Dest>
00217   void _solveWithGuess(const Rhs& b, Dest& x) const
  {    
    for(int j=0; j<b.cols(); ++j)
    {
      m_iterations = Base::m_maxIterations;
      m_error = Base::m_tolerance;
      
      typename Dest::ColXpr xj(x,j);
      internal::bicgstab(*mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_error);
    }
    
    m_isInitialized = true;
    m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
  }

  /** \internal */
  template<typename Rhs,typename Dest>
00234   void _solve(const Rhs& b, Dest& x) const
  {
    x.setOnes();
    _solveWithGuess(b,x);
  }

protected:

};


namespace internal {

  template<typename _MatrixType, typename _Preconditioner, typename Rhs>
00248 struct solve_retval<BiCGSTAB<_MatrixType, _Preconditioner>, Rhs>
  : solve_retval_base<BiCGSTAB<_MatrixType, _Preconditioner>, Rhs>
{
  typedef BiCGSTAB<_MatrixType, _Preconditioner> Dec;
  EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)

  template<typename Dest> void evalTo(Dest& dst) const
  {
    dec()._solve(rhs(),dst);
  }
};

}

#endif // EIGEN_BICGSTAB_H

Generated by  Doxygen 1.6.0   Back to index