13 #include "./InternalHeaderCheck.h"
17 template<
typename Decomposition,
typename RhsType,
typename StorageKind>
class SolveImpl;
34 template<
typename Decomposition,
typename RhsType,
typename StorageKind>
struct solve_traits;
36 template<
typename Decomposition,
typename RhsType>
37 struct solve_traits<Decomposition,RhsType,Dense>
39 typedef typename make_proper_matrix_type<
typename RhsType::Scalar,
40 Decomposition::ColsAtCompileTime,
41 RhsType::ColsAtCompileTime,
42 RhsType::PlainObject::Options,
43 Decomposition::MaxColsAtCompileTime,
44 RhsType::MaxColsAtCompileTime>::type PlainObject;
47 template<
typename Decomposition,
typename RhsType>
48 struct traits<Solve<Decomposition, RhsType> >
49 : traits<typename solve_traits<Decomposition,RhsType,typename internal::traits<RhsType>::StorageKind>::PlainObject>
51 typedef typename solve_traits<Decomposition,RhsType,typename internal::traits<RhsType>::StorageKind>::PlainObject PlainObject;
52 typedef typename promote_index_type<typename Decomposition::StorageIndex, typename RhsType::StorageIndex>::type StorageIndex;
53 typedef traits<PlainObject> BaseTraits;
63 template<
typename Decomposition,
typename RhsType>
64 class Solve :
public SolveImpl<Decomposition,RhsType,typename internal::traits<RhsType>::StorageKind>
67 typedef typename internal::traits<Solve>::PlainObject PlainObject;
68 typedef typename internal::traits<Solve>::StorageIndex StorageIndex;
70 Solve(
const Decomposition &dec,
const RhsType &rhs)
71 : m_dec(dec), m_rhs(rhs)
74 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
Index rows()
const EIGEN_NOEXCEPT {
return m_dec.cols(); }
75 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
Index cols()
const EIGEN_NOEXCEPT {
return m_rhs.cols(); }
77 EIGEN_DEVICE_FUNC
const Decomposition& dec()
const {
return m_dec; }
78 EIGEN_DEVICE_FUNC
const RhsType& rhs()
const {
return m_rhs; }
81 const Decomposition &m_dec;
82 const typename internal::ref_selector<RhsType>::type m_rhs;
87 template<
typename Decomposition,
typename RhsType>
88 class SolveImpl<Decomposition,RhsType,
Dense>
89 :
public MatrixBase<Solve<Decomposition,RhsType> >
96 EIGEN_DENSE_PUBLIC_INTERFACE(Derived)
101 Scalar coeff(
Index i)
const;
105 template<
typename Decomposition,
typename RhsType,
typename StorageKind>
106 class SolveImpl :
public internal::generic_xpr_base<Solve<Decomposition,RhsType>, MatrixXpr, StorageKind>::type
109 typedef typename internal::generic_xpr_base<Solve<Decomposition,RhsType>, MatrixXpr, StorageKind>::type Base;
115 template<
typename Decomposition,
typename RhsType>
116 struct evaluator<Solve<Decomposition,RhsType> >
117 :
public evaluator<typename Solve<Decomposition,RhsType>::PlainObject>
119 typedef Solve<Decomposition,RhsType> SolveType;
120 typedef typename SolveType::PlainObject PlainObject;
121 typedef evaluator<PlainObject> Base;
125 EIGEN_DEVICE_FUNC
explicit evaluator(
const SolveType& solve)
126 : m_result(solve.rows(), solve.cols())
128 internal::construct_at<Base>(
this, m_result);
129 solve.dec()._solve_impl(solve.rhs(), m_result);
133 PlainObject m_result;
138 template<
typename DstXprType,
typename DecType,
typename RhsType,
typename Scalar>
139 struct Assignment<DstXprType, Solve<DecType,RhsType>, internal::assign_op<Scalar,Scalar>, Dense2Dense>
141 typedef Solve<DecType,RhsType> SrcXprType;
142 static void run(DstXprType &dst,
const SrcXprType &src,
const internal::assign_op<Scalar,Scalar> &)
144 Index dstRows = src.rows();
145 Index dstCols = src.cols();
146 if((dst.rows()!=dstRows) || (dst.cols()!=dstCols))
147 dst.resize(dstRows, dstCols);
149 src.dec()._solve_impl(src.rhs(), dst);
154 template<
typename DstXprType,
typename DecType,
typename RhsType,
typename Scalar>
155 struct Assignment<DstXprType, Solve<Transpose<const DecType>,RhsType>, internal::assign_op<Scalar,Scalar>, Dense2Dense>
157 typedef Solve<Transpose<const DecType>,RhsType> SrcXprType;
158 static void run(DstXprType &dst,
const SrcXprType &src,
const internal::assign_op<Scalar,Scalar> &)
160 Index dstRows = src.rows();
161 Index dstCols = src.cols();
162 if((dst.rows()!=dstRows) || (dst.cols()!=dstCols))
163 dst.resize(dstRows, dstCols);
165 src.dec().nestedExpression().template _solve_impl_transposed<false>(src.rhs(), dst);
170 template<
typename DstXprType,
typename DecType,
typename RhsType,
typename Scalar>
171 struct Assignment<DstXprType, Solve<CwiseUnaryOp<internal::scalar_conjugate_op<typename DecType::Scalar>, const Transpose<const DecType> >,RhsType>,
172 internal::assign_op<Scalar,Scalar>, Dense2Dense>
174 typedef Solve<CwiseUnaryOp<internal::scalar_conjugate_op<typename DecType::Scalar>,
const Transpose<const DecType> >,RhsType> SrcXprType;
175 static void run(DstXprType &dst,
const SrcXprType &src,
const internal::assign_op<Scalar,Scalar> &)
177 Index dstRows = src.rows();
178 Index dstCols = src.cols();
179 if((dst.rows()!=dstRows) || (dst.cols()!=dstCols))
180 dst.resize(dstRows, dstCols);
182 src.dec().nestedExpression().nestedExpression().template _solve_impl_transposed<true>(src.rhs(), dst);
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
Pseudo expression representing a solving operation.
Definition: Solve.h:65
const unsigned int EvalBeforeNestingBit
Definition: Constants.h:72
const unsigned int RowMajorBit
Definition: Constants.h:68
Namespace containing all symbols from the Eigen library.
Definition: Core:139
const int HugeCost
Definition: Constants.h:46
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
Definition: Constants.h:509