10 #ifndef EIGEN_TRIANGULAR_SOLVER_VECTOR_H
11 #define EIGEN_TRIANGULAR_SOLVER_VECTOR_H
13 #include "../InternalHeaderCheck.h"
19 template<
typename LhsScalar,
typename RhsScalar,
typename Index,
int Mode,
bool Conjugate,
int StorageOrder>
20 struct triangular_solve_vector<LhsScalar, RhsScalar,
Index,
OnTheRight, Mode, Conjugate, StorageOrder>
22 static void run(
Index size,
const LhsScalar* _lhs,
Index lhsStride, RhsScalar* rhs)
27 >::run(
size, _lhs, lhsStride, rhs);
32 template<
typename LhsScalar,
typename RhsScalar,
typename Index,
int Mode,
bool Conjugate>
38 static void run(
Index size,
const LhsScalar* _lhs,
Index lhsStride, RhsScalar* rhs)
40 typedef Map<const Matrix<LhsScalar,Dynamic,Dynamic,RowMajor>, 0, OuterStride<> > LhsMap;
41 const LhsMap lhs(_lhs,
size,
size,OuterStride<>(lhsStride));
43 typedef const_blas_data_mapper<LhsScalar,Index,RowMajor> LhsMapper;
44 typedef const_blas_data_mapper<RhsScalar,Index,ColMajor> RhsMapper;
48 const CwiseUnaryOp<typename internal::scalar_conjugate_op<LhsScalar>,LhsMap>,
49 const LhsMap&> cjLhs(lhs);
52 IsLower ? pi<size : pi>0;
53 IsLower ? pi+=PanelWidth : pi-=PanelWidth)
63 Index startRow = IsLower ? pi : pi-actualPanelWidth;
64 Index startCol = IsLower ? 0 : pi;
66 general_matrix_vector_product<Index,LhsScalar,LhsMapper,RowMajor,Conjugate,RhsScalar,RhsMapper,false>::run(
68 LhsMapper(&lhs.coeffRef(startRow,startCol), lhsStride),
69 RhsMapper(rhs + startCol, 1),
74 for(
Index k=0; k<actualPanelWidth; ++k)
76 Index i = IsLower ? pi+k : pi-k-1;
77 Index s = IsLower ? pi :
i+1;
79 rhs[
i] -= (cjLhs.row(
i).segment(s,k).transpose().cwiseProduct(Map<
const Matrix<RhsScalar,Dynamic,1> >(rhs+s,k))).sum();
89 template<
typename LhsScalar,
typename RhsScalar,
typename Index,
int Mode,
bool Conjugate>
95 static void run(
Index size,
const LhsScalar* _lhs,
Index lhsStride, RhsScalar* rhs)
97 typedef Map<const Matrix<LhsScalar,Dynamic,Dynamic,ColMajor>, 0, OuterStride<> > LhsMap;
98 const LhsMap lhs(_lhs,
size,
size,OuterStride<>(lhsStride));
99 typedef const_blas_data_mapper<LhsScalar,Index,ColMajor> LhsMapper;
100 typedef const_blas_data_mapper<RhsScalar,Index,ColMajor> RhsMapper;
101 std::conditional_t<Conjugate,
102 const CwiseUnaryOp<typename internal::scalar_conjugate_op<LhsScalar>,LhsMap>,
108 IsLower ? pi<size : pi>0;
109 IsLower ? pi+=PanelWidth : pi-=PanelWidth)
112 Index startBlock = IsLower ? pi : pi-actualPanelWidth;
113 Index endBlock = IsLower ? pi + actualPanelWidth : 0;
115 for(
Index k=0; k<actualPanelWidth; ++k)
117 Index i = IsLower ? pi+k : pi-k-1;
121 rhs[
i] /= cjLhs.coeff(
i,
i);
123 Index r = actualPanelWidth - k - 1;
126 Map<Matrix<RhsScalar,Dynamic,1> >(rhs+s,r) -= rhs[
i] * cjLhs.col(
i).segment(s,r);
129 Index r = IsLower ?
size - endBlock : startBlock;
135 general_matrix_vector_product<Index,LhsScalar,LhsMapper,ColMajor,Conjugate,RhsScalar,RhsMapper,false>::run(
137 LhsMapper(&lhs.coeffRef(endBlock,startBlock), lhsStride),
138 RhsMapper(rhs+startBlock, 1),
139 rhs+endBlock, 1, RhsScalar(-1));
#define EIGEN_TUNE_TRIANGULAR_PANEL_WIDTH
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
bool is_identically_zero(const Scalar &s)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.