Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/eigen3/unsupported/Eigen/src/IterativeSolvers/IDRS.h was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2020 Chris Schoutrop <c.e.m.schoutrop@tue.nl>
0005 // Copyright (C) 2020 Jens Wehner <j.wehner@esciencecenter.nl>
0006 // Copyright (C) 2020 Jan van Dijk <j.v.dijk@tue.nl>
0007 //
0008 // This Source Code Form is subject to the terms of the Mozilla
0009 // Public License v. 2.0. If a copy of the MPL was not distributed
0010 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0011 
0012 
0013 #ifndef EIGEN_IDRS_H
0014 #define EIGEN_IDRS_H
0015 
0016 namespace Eigen
0017 {
0018 
0019     namespace internal
0020     {
0021         /**     \internal Low-level Induced Dimension Reduction algoritm
0022                 \param A The matrix A
0023                 \param b The right hand side vector b
0024                 \param x On input and initial solution, on output the computed solution.
0025                 \param precond A preconditioner being able to efficiently solve for an
0026                           approximation of Ax=b (regardless of b)
0027                 \param iter On input the max number of iteration, on output the number of performed iterations.
0028                 \param relres On input the tolerance error, on output an estimation of the relative error.
0029                 \param S On input Number of the dimension of the shadow space.
0030                 \param smoothing switches residual smoothing on.
0031                 \param angle small omega lead to faster convergence at the expense of numerical stability
0032                 \param replacement switches on a residual replacement strategy to increase accuracy of residual at the expense of more Mat*vec products
0033                 \return false in the case of numerical issue, for example a break down of IDRS.
0034         */
0035         template<typename Vector, typename RealScalar>
0036         typename Vector::Scalar omega(const Vector& t, const Vector& s, RealScalar angle)
0037         {
0038             using numext::abs;
0039             typedef typename Vector::Scalar Scalar;
0040             const RealScalar ns = s.norm();
0041             const RealScalar nt = t.norm();
0042             const Scalar ts = t.dot(s);
0043             const RealScalar rho = abs(ts / (nt * ns));
0044 
0045             if (rho < angle) {
0046                 if (ts == Scalar(0)) {
0047                     return Scalar(0);
0048                 }
0049                 // Original relation for om is given by
0050                 // om = om * angle / rho;
0051                 // To alleviate potential (near) division by zero this can be rewritten as
0052                 // om = angle * (ns / nt) * (ts / abs(ts)) = angle * (ns / nt) * sgn(ts)
0053                 return angle * (ns / nt) * (ts / abs(ts));
0054             }
0055             return ts / (nt * nt);
0056         }
0057 
0058         template <typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
0059         bool idrs(const MatrixType& A, const Rhs& b, Dest& x, const Preconditioner& precond,
0060             Index& iter,
0061             typename Dest::RealScalar& relres, Index S, bool smoothing, typename Dest::RealScalar angle, bool replacement)
0062         {
0063             typedef typename Dest::RealScalar RealScalar;
0064             typedef typename Dest::Scalar Scalar;
0065             typedef Matrix<Scalar, Dynamic, 1> VectorType;
0066             typedef Matrix<Scalar, Dynamic, Dynamic, ColMajor> DenseMatrixType;
0067             const Index N = b.size();
0068             S = S < x.rows() ? S : x.rows();
0069             const RealScalar tol = relres;
0070             const Index maxit = iter;
0071 
0072             Index replacements = 0;
0073             bool trueres = false;
0074 
0075             FullPivLU<DenseMatrixType> lu_solver;
0076 
0077             DenseMatrixType P;
0078             {
0079                 HouseholderQR<DenseMatrixType> qr(DenseMatrixType::Random(N, S));
0080                 P = (qr.householderQ() * DenseMatrixType::Identity(N, S));
0081             }
0082 
0083             const RealScalar normb = b.norm();
0084 
0085             if (internal::isApprox(normb, RealScalar(0)))
0086             {
0087                 //Solution is the zero vector
0088                 x.setZero();
0089                 iter = 0;
0090                 relres = 0;
0091                 return true;
0092             }
0093              // from http://homepage.tudelft.nl/1w5b5/IDRS/manual.pdf
0094              // A peak in the residual is considered dangerously high if‖ri‖/‖b‖> C(tol/epsilon).
0095              // With epsilon the
0096              // relative machine precision. The factor tol/epsilon corresponds to the size of a
0097              // finite precision number that is so large that the absolute round-off error in
0098              // this number, when propagated through the process, makes it impossible to
0099              // achieve the required accuracy.The factor C accounts for the accumulation of
0100              // round-off errors. This parameter has beenset to 10−3.
0101              // mp is epsilon/C
0102              // 10^3 * eps is very conservative, so normally no residual replacements will take place. 
0103              // It only happens if things go very wrong. Too many restarts may ruin the convergence.
0104             const RealScalar mp = RealScalar(1e3) * NumTraits<Scalar>::epsilon();
0105 
0106 
0107 
0108             //Compute initial residual
0109             const RealScalar tolb = tol * normb; //Relative tolerance
0110             VectorType r = b - A * x;
0111 
0112             VectorType x_s, r_s;
0113 
0114             if (smoothing)
0115             {
0116                 x_s = x;
0117                 r_s = r;
0118             }
0119 
0120             RealScalar normr = r.norm();
0121 
0122             if (normr <= tolb)
0123             {
0124                 //Initial guess is a good enough solution
0125                 iter = 0;
0126                 relres = normr / normb;
0127                 return true;
0128             }
0129 
0130             DenseMatrixType G = DenseMatrixType::Zero(N, S);
0131             DenseMatrixType U = DenseMatrixType::Zero(N, S);
0132             DenseMatrixType M = DenseMatrixType::Identity(S, S);
0133             VectorType t(N), v(N);
0134             Scalar om = 1.;
0135 
0136             //Main iteration loop, guild G-spaces:
0137             iter = 0;
0138 
0139             while (normr > tolb && iter < maxit)
0140             {
0141                 //New right hand size for small system:
0142                 VectorType f = (r.adjoint() * P).adjoint();
0143 
0144                 for (Index k = 0; k < S; ++k)
0145                 {
0146                     //Solve small system and make v orthogonal to P:
0147                     //c = M(k:s,k:s)\f(k:s);
0148                     lu_solver.compute(M.block(k , k , S -k, S - k ));
0149                     VectorType c = lu_solver.solve(f.segment(k , S - k ));
0150                     //v = r - G(:,k:s)*c;
0151                     v = r - G.rightCols(S - k ) * c;
0152                     //Preconditioning
0153                     v = precond.solve(v);
0154 
0155                     //Compute new U(:,k) and G(:,k), G(:,k) is in space G_j
0156                     U.col(k) = U.rightCols(S - k ) * c + om * v;
0157                     G.col(k) = A * U.col(k );
0158 
0159                     //Bi-Orthogonalise the new basis vectors:
0160                     for (Index i = 0; i < k-1 ; ++i)
0161                     {
0162                         //alpha =  ( P(:,i)'*G(:,k) )/M(i,i);
0163                         Scalar alpha = P.col(i ).dot(G.col(k )) / M(i, i );
0164                         G.col(k ) = G.col(k ) - alpha * G.col(i );
0165                         U.col(k ) = U.col(k ) - alpha * U.col(i );
0166                     }
0167 
0168                     //New column of M = P'*G  (first k-1 entries are zero)
0169                     //M(k:s,k) = (G(:,k)'*P(:,k:s))';
0170                     M.block(k , k , S - k , 1) = (G.col(k ).adjoint() * P.rightCols(S - k )).adjoint();
0171 
0172                     if (internal::isApprox(M(k,k), Scalar(0)))
0173                     {
0174                         return false;
0175                     }
0176 
0177                     //Make r orthogonal to q_i, i = 0..k-1
0178                     Scalar beta = f(k ) / M(k , k );
0179                     r = r - beta * G.col(k );
0180                     x = x + beta * U.col(k );
0181                     normr = r.norm();
0182 
0183                     if (replacement && normr > tolb / mp)
0184                     {
0185                         trueres = true;
0186                     }
0187 
0188                     //Smoothing:
0189                     if (smoothing)
0190                     {
0191                         t = r_s - r;
0192                         //gamma is a Scalar, but the conversion is not allowed
0193                         Scalar gamma = t.dot(r_s) / t.norm();
0194                         r_s = r_s - gamma * t;
0195                         x_s = x_s - gamma * (x_s - x);
0196                         normr = r_s.norm();
0197                     }
0198 
0199                     if (normr < tolb || iter == maxit)
0200                     {
0201                         break;
0202                     }
0203 
0204                     //New f = P'*r (first k  components are zero)
0205                     if (k < S-1)
0206                     {
0207                         f.segment(k + 1, S - (k + 1) ) = f.segment(k + 1 , S - (k + 1)) - beta * M.block(k + 1 , k , S - (k + 1), 1);
0208                     }
0209                 }//end for
0210 
0211                 if (normr < tolb || iter == maxit)
0212                 {
0213                     break;
0214                 }
0215 
0216                 //Now we have sufficient vectors in G_j to compute residual in G_j+1
0217                 //Note: r is already perpendicular to P so v = r
0218                 //Preconditioning
0219                 v = r;
0220                 v = precond.solve(v);
0221 
0222                 //Matrix-vector multiplication:
0223                 t = A * v;
0224 
0225                 //Computation of a new omega
0226                 om = internal::omega(t, r, angle);
0227 
0228                 if (om == RealScalar(0.0))
0229                 {
0230                     return false;
0231                 }
0232 
0233                 r = r - om * t;
0234                 x = x + om * v;
0235                 normr = r.norm();
0236 
0237                 if (replacement && normr > tolb / mp)
0238                 {
0239                     trueres = true;
0240                 }
0241 
0242                 //Residual replacement?
0243                 if (trueres && normr < normb)
0244                 {
0245                     r = b - A * x;
0246                     trueres = false;
0247                     replacements++;
0248                 }
0249 
0250                 //Smoothing:
0251                 if (smoothing)
0252                 {
0253                     t = r_s - r;
0254                     Scalar gamma = t.dot(r_s) /t.norm();
0255                     r_s = r_s - gamma * t;
0256                     x_s = x_s - gamma * (x_s - x);
0257                     normr = r_s.norm();
0258                 }
0259 
0260                 iter++;
0261 
0262             }//end while
0263 
0264             if (smoothing)
0265             {
0266                 x = x_s;
0267             }
0268             relres=normr/normb;
0269             return true;
0270         }
0271 
0272     }  // namespace internal
0273 
0274     template <typename _MatrixType, typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
0275     class IDRS;
0276 
0277     namespace internal
0278     {
0279 
0280         template <typename _MatrixType, typename _Preconditioner>
0281         struct traits<Eigen::IDRS<_MatrixType, _Preconditioner> >
0282         {
0283             typedef _MatrixType MatrixType;
0284             typedef _Preconditioner Preconditioner;
0285         };
0286 
0287     }  // namespace internal
0288 
0289 
0290 /** \ingroup IterativeLinearSolvers_Module
0291   * \brief The Induced Dimension Reduction method (IDR(s)) is a short-recurrences Krylov method for sparse square problems.
0292   *
0293   * This class allows to solve for A.x = b sparse linear problems. The vectors x and b can be either dense or sparse.
0294   * he Induced Dimension Reduction method, IDR(), is a robust and efficient short-recurrence Krylov subspace method for
0295   * solving large nonsymmetric systems of linear equations.
0296   *
0297   * For indefinite systems IDR(S) outperforms both BiCGStab and BiCGStab(L). Additionally, IDR(S) can handle matrices
0298   * with complex eigenvalues more efficiently than BiCGStab.
0299   *
0300   * Many problems that do not converge for BiCGSTAB converge for IDR(s) (for larger values of s). And if both methods 
0301   * converge the convergence for IDR(s) is typically much faster for difficult systems (for example indefinite problems). 
0302   *
0303   * IDR(s) is a limited memory finite termination method. In exact arithmetic it converges in at most N+N/s iterations,
0304   * with N the system size.  It uses a fixed number of 4+3s vector. In comparison, BiCGSTAB terminates in 2N iterations 
0305   * and uses 7 vectors. GMRES terminates in at most N iterations, and uses I+3 vectors, with I the number of iterations. 
0306   * Restarting GMRES limits the memory consumption, but destroys the finite termination property.
0307   *
0308   * \tparam _MatrixType the type of the sparse matrix A, can be a dense or a sparse matrix.
0309   * \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
0310   *
0311   * \implsparsesolverconcept
0312   *
0313   * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
0314   * and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
0315   * and NumTraits<Scalar>::epsilon() for the tolerance.
0316   *
0317   * The tolerance corresponds to the relative residual error: |Ax-b|/|b|
0318   *
0319   * \b Performance: when using sparse matrices, best performance is achied for a row-major sparse matrix format.
0320   * Moreover, in this case multi-threading can be exploited if the user code is compiled with OpenMP enabled.
0321   * See \ref TopicMultiThreading for details.
0322   *
0323   * By default the iterations start with x=0 as an initial guess of the solution.
0324   * One can control the start using the solveWithGuess() method.
0325   *
0326   * IDR(s) can also be used in a matrix-free context, see the following \link MatrixfreeSolverExample example \endlink.
0327   *
0328   * \sa class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
0329   */
0330     template <typename _MatrixType, typename _Preconditioner>
0331     class IDRS : public IterativeSolverBase<IDRS<_MatrixType, _Preconditioner> >
0332     {
0333 
0334         public:
0335             typedef _MatrixType MatrixType;
0336             typedef typename MatrixType::Scalar Scalar;
0337             typedef typename MatrixType::RealScalar RealScalar;
0338             typedef _Preconditioner Preconditioner;
0339 
0340         private:
0341             typedef IterativeSolverBase<IDRS> Base;
0342             using Base::m_error;
0343             using Base::m_info;
0344             using Base::m_isInitialized;
0345             using Base::m_iterations;
0346             using Base::matrix;
0347             Index m_S;
0348             bool m_smoothing;
0349             RealScalar m_angle;
0350             bool m_residual;
0351 
0352         public:
0353             /** Default constructor. */
0354             IDRS(): m_S(4), m_smoothing(false), m_angle(RealScalar(0.7)), m_residual(false) {}
0355 
0356             /**     Initialize the solver with matrix \a A for further \c Ax=b solving.
0357 
0358                     This constructor is a shortcut for the default constructor followed
0359                     by a call to compute().
0360 
0361                     \warning this class stores a reference to the matrix A as well as some
0362                     precomputed values that depend on it. Therefore, if \a A is changed
0363                     this class becomes invalid. Call compute() to update it with the new
0364                     matrix A, or modify a copy of A.
0365             */
0366             template <typename MatrixDerived>
0367             explicit IDRS(const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_S(4), m_smoothing(false),
0368                                                                m_angle(RealScalar(0.7)), m_residual(false) {}
0369 
0370 
0371             /** \internal */
0372             /**     Loops over the number of columns of b and does the following:
0373                             1. sets the tolerence and maxIterations
0374                             2. Calls the function that has the core solver routine
0375             */
0376             template <typename Rhs, typename Dest>
0377             void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const
0378             {
0379                 m_iterations = Base::maxIterations();
0380                 m_error = Base::m_tolerance;
0381 
0382                 bool ret = internal::idrs(matrix(), b, x, Base::m_preconditioner, m_iterations, m_error, m_S,m_smoothing,m_angle,m_residual);
0383 
0384                 m_info = (!ret) ? NumericalIssue : m_error <= Base::m_tolerance ? Success : NoConvergence;
0385             }
0386 
0387             /** Sets the parameter S, indicating the dimension of the shadow space. Default is 4*/
0388             void setS(Index S)
0389             {
0390                 if (S < 1)
0391                 {
0392                     S = 4;
0393                 }
0394 
0395                 m_S = S;
0396             }
0397 
0398             /** Switches off and on smoothing.
0399             Residual smoothing results in monotonically decreasing residual norms at
0400             the expense of two extra vectors of storage and a few extra vector
0401             operations. Although monotonic decrease of the residual norms is a
0402             desirable property, the rate of convergence of the unsmoothed process and
0403             the smoothed process is basically the same. Default is off */
0404             void setSmoothing(bool smoothing)
0405             {
0406                 m_smoothing=smoothing;
0407             }
0408 
0409             /** The angle must be a real scalar. In IDR(s), a value for the
0410             iteration parameter omega must be chosen in every s+1th step. The most
0411             natural choice is to select a value to minimize the norm of the next residual.
0412             This corresponds to the parameter omega = 0. In practice, this may lead to
0413             values of omega that are so small that the other iteration parameters
0414             cannot be computed with sufficient accuracy. In such cases it is better to
0415             increase the value of omega sufficiently such that a compromise is reached
0416             between accurate computations and reduction of the residual norm. The
0417             parameter angle =0.7 (”maintaining the convergence strategy”)
0418             results in such a compromise. */
0419             void setAngle(RealScalar angle)
0420             {
0421                 m_angle=angle;
0422             }
0423 
0424             /** The parameter replace is a logical that determines whether a
0425             residual replacement strategy is employed to increase the accuracy of the
0426             solution. */
0427             void setResidualUpdate(bool update)
0428             {
0429                 m_residual=update;
0430             }
0431 
0432     };
0433 
0434 }  // namespace Eigen
0435 
0436 #endif /* EIGEN_IDRS_H */