xref: /aosp_15_r20/external/eigen/Eigen/src/IterativeLinearSolvers/SolveWithGuess.h (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2014 Gael Guennebaud <[email protected]>
5*bf2c3715SXin Li //
6*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
7*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
8*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9*bf2c3715SXin Li 
10*bf2c3715SXin Li #ifndef EIGEN_SOLVEWITHGUESS_H
11*bf2c3715SXin Li #define EIGEN_SOLVEWITHGUESS_H
12*bf2c3715SXin Li 
13*bf2c3715SXin Li namespace Eigen {
14*bf2c3715SXin Li 
15*bf2c3715SXin Li template<typename Decomposition, typename RhsType, typename GuessType> class SolveWithGuess;
16*bf2c3715SXin Li 
17*bf2c3715SXin Li /** \class SolveWithGuess
18*bf2c3715SXin Li   * \ingroup IterativeLinearSolvers_Module
19*bf2c3715SXin Li   *
20*bf2c3715SXin Li   * \brief Pseudo expression representing a solving operation
21*bf2c3715SXin Li   *
22*bf2c3715SXin Li   * \tparam Decomposition the type of the matrix or decomposion object
23*bf2c3715SXin Li   * \tparam Rhstype the type of the right-hand side
24*bf2c3715SXin Li   *
25*bf2c3715SXin Li   * This class represents an expression of A.solve(B)
26*bf2c3715SXin Li   * and most of the time this is the only way it is used.
27*bf2c3715SXin Li   *
28*bf2c3715SXin Li   */
29*bf2c3715SXin Li namespace internal {
30*bf2c3715SXin Li 
31*bf2c3715SXin Li 
32*bf2c3715SXin Li template<typename Decomposition, typename RhsType, typename GuessType>
33*bf2c3715SXin Li struct traits<SolveWithGuess<Decomposition, RhsType, GuessType> >
34*bf2c3715SXin Li   : traits<Solve<Decomposition,RhsType> >
35*bf2c3715SXin Li {};
36*bf2c3715SXin Li 
37*bf2c3715SXin Li }
38*bf2c3715SXin Li 
39*bf2c3715SXin Li 
40*bf2c3715SXin Li template<typename Decomposition, typename RhsType, typename GuessType>
41*bf2c3715SXin Li class SolveWithGuess : public internal::generic_xpr_base<SolveWithGuess<Decomposition,RhsType,GuessType>, MatrixXpr, typename internal::traits<RhsType>::StorageKind>::type
42*bf2c3715SXin Li {
43*bf2c3715SXin Li public:
44*bf2c3715SXin Li   typedef typename internal::traits<SolveWithGuess>::Scalar Scalar;
45*bf2c3715SXin Li   typedef typename internal::traits<SolveWithGuess>::PlainObject PlainObject;
46*bf2c3715SXin Li   typedef typename internal::generic_xpr_base<SolveWithGuess<Decomposition,RhsType,GuessType>, MatrixXpr, typename internal::traits<RhsType>::StorageKind>::type Base;
47*bf2c3715SXin Li   typedef typename internal::ref_selector<SolveWithGuess>::type Nested;
48*bf2c3715SXin Li 
49*bf2c3715SXin Li   SolveWithGuess(const Decomposition &dec, const RhsType &rhs, const GuessType &guess)
50*bf2c3715SXin Li     : m_dec(dec), m_rhs(rhs), m_guess(guess)
51*bf2c3715SXin Li   {}
52*bf2c3715SXin Li 
53*bf2c3715SXin Li   EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
54*bf2c3715SXin Li   Index rows() const EIGEN_NOEXCEPT { return m_dec.cols(); }
55*bf2c3715SXin Li   EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
56*bf2c3715SXin Li   Index cols() const EIGEN_NOEXCEPT { return m_rhs.cols(); }
57*bf2c3715SXin Li 
58*bf2c3715SXin Li   EIGEN_DEVICE_FUNC const Decomposition& dec()   const { return m_dec; }
59*bf2c3715SXin Li   EIGEN_DEVICE_FUNC const RhsType&       rhs()   const { return m_rhs; }
60*bf2c3715SXin Li   EIGEN_DEVICE_FUNC const GuessType&     guess() const { return m_guess; }
61*bf2c3715SXin Li 
62*bf2c3715SXin Li protected:
63*bf2c3715SXin Li   const Decomposition &m_dec;
64*bf2c3715SXin Li   const RhsType       &m_rhs;
65*bf2c3715SXin Li   const GuessType     &m_guess;
66*bf2c3715SXin Li 
67*bf2c3715SXin Li private:
68*bf2c3715SXin Li   Scalar coeff(Index row, Index col) const;
69*bf2c3715SXin Li   Scalar coeff(Index i) const;
70*bf2c3715SXin Li };
71*bf2c3715SXin Li 
72*bf2c3715SXin Li namespace internal {
73*bf2c3715SXin Li 
74*bf2c3715SXin Li // Evaluator of SolveWithGuess -> eval into a temporary
75*bf2c3715SXin Li template<typename Decomposition, typename RhsType, typename GuessType>
76*bf2c3715SXin Li struct evaluator<SolveWithGuess<Decomposition,RhsType, GuessType> >
77*bf2c3715SXin Li   : public evaluator<typename SolveWithGuess<Decomposition,RhsType,GuessType>::PlainObject>
78*bf2c3715SXin Li {
79*bf2c3715SXin Li   typedef SolveWithGuess<Decomposition,RhsType,GuessType> SolveType;
80*bf2c3715SXin Li   typedef typename SolveType::PlainObject PlainObject;
81*bf2c3715SXin Li   typedef evaluator<PlainObject> Base;
82*bf2c3715SXin Li 
83*bf2c3715SXin Li   evaluator(const SolveType& solve)
84*bf2c3715SXin Li     : m_result(solve.rows(), solve.cols())
85*bf2c3715SXin Li   {
86*bf2c3715SXin Li     ::new (static_cast<Base*>(this)) Base(m_result);
87*bf2c3715SXin Li     m_result = solve.guess();
88*bf2c3715SXin Li     solve.dec()._solve_with_guess_impl(solve.rhs(), m_result);
89*bf2c3715SXin Li   }
90*bf2c3715SXin Li 
91*bf2c3715SXin Li protected:
92*bf2c3715SXin Li   PlainObject m_result;
93*bf2c3715SXin Li };
94*bf2c3715SXin Li 
95*bf2c3715SXin Li // Specialization for "dst = dec.solveWithGuess(rhs)"
96*bf2c3715SXin Li // NOTE we need to specialize it for Dense2Dense to avoid ambiguous specialization error and a Sparse2Sparse specialization must exist somewhere
97*bf2c3715SXin Li template<typename DstXprType, typename DecType, typename RhsType, typename GuessType, typename Scalar>
98*bf2c3715SXin Li struct Assignment<DstXprType, SolveWithGuess<DecType,RhsType,GuessType>, internal::assign_op<Scalar,Scalar>, Dense2Dense>
99*bf2c3715SXin Li {
100*bf2c3715SXin Li   typedef SolveWithGuess<DecType,RhsType,GuessType> SrcXprType;
101*bf2c3715SXin Li   static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<Scalar,Scalar> &)
102*bf2c3715SXin Li   {
103*bf2c3715SXin Li     Index dstRows = src.rows();
104*bf2c3715SXin Li     Index dstCols = src.cols();
105*bf2c3715SXin Li     if((dst.rows()!=dstRows) || (dst.cols()!=dstCols))
106*bf2c3715SXin Li       dst.resize(dstRows, dstCols);
107*bf2c3715SXin Li 
108*bf2c3715SXin Li     dst = src.guess();
109*bf2c3715SXin Li     src.dec()._solve_with_guess_impl(src.rhs(), dst/*, src.guess()*/);
110*bf2c3715SXin Li   }
111*bf2c3715SXin Li };
112*bf2c3715SXin Li 
113*bf2c3715SXin Li } // end namespace internal
114*bf2c3715SXin Li 
115*bf2c3715SXin Li } // end namespace Eigen
116*bf2c3715SXin Li 
117*bf2c3715SXin Li #endif // EIGEN_SOLVEWITHGUESS_H
118