11#ifndef EIGEN_TRIANGULAR_SOLVER_MATRIX_H
12#define EIGEN_TRIANGULAR_SOLVER_MATRIX_H
15#include "../InternalHeaderCheck.h"
21template <
typename Scalar,
typename Index,
int Mode,
bool Conjugate,
int TriStorageOrder,
int OtherInnerStride,
26 static void kernel(
Index size,
Index otherSize,
const Scalar* _tri,
Index triStride, Scalar* _other,
Index otherIncr,
30template <
typename Scalar,
typename Index,
int Mode,
bool Conjugate,
int TriStorageOrder,
int OtherInnerStride,
35 static void kernel(
Index size,
Index otherSize,
const Scalar* _tri,
Index triStride, Scalar* _other,
Index otherIncr,
39template <
typename Scalar,
typename Index,
int Mode,
bool Conjugate,
int TriStorageOrder,
int OtherInnerStride,
41EIGEN_STRONG_INLINE
void trsmKernelL<Scalar,
Index, Mode, Conjugate, TriStorageOrder, OtherInnerStride,
42 Specialized>::kernel(
Index size,
Index otherSize,
const Scalar* _tri,
43 Index triStride, Scalar* _other,
Index otherIncr,
45 typedef const_blas_data_mapper<Scalar, Index, TriStorageOrder> TriMapper;
46 typedef blas_data_mapper<Scalar, Index, ColMajor, Unaligned, OtherInnerStride> OtherMapper;
47 TriMapper tri(_tri, triStride);
48 OtherMapper other(_other, otherStride, otherIncr);
51 conj_if<Conjugate>
conj;
54 for (
Index k = 0; k < size; ++k) {
56 Index i = IsLower ? k : -k - 1;
57 Index rs = size - k - 1;
58 Index s = TriStorageOrder ==
RowMajor ? (IsLower ? 0 : i + 1) : IsLower ? i + 1 : i - rs;
60 Scalar a = (Mode &
UnitDiag) ? Scalar(1) : Scalar(Scalar(1) /
conj(tri(i, i)));
61 for (
Index j = 0; j < otherSize; ++j) {
64 const Scalar* l = &tri(i, s);
65 typename OtherMapper::LinearMapper r = other.getLinearMapper(s, j);
66 for (
Index i3 = 0; i3 < k; ++i3) b +=
conj(l[i3]) * r(i3);
68 other(i, j) = (other(i, j) - b) * a;
70 Scalar& otherij = other(i, j);
73 typename OtherMapper::LinearMapper r = other.getLinearMapper(s, j);
74 typename TriMapper::LinearMapper l = tri.getLinearMapper(s, i);
75 for (
Index i3 = 0; i3 < rs; ++i3) r(i3) -= b *
conj(l(i3));
81template <
typename Scalar,
typename Index,
int Mode,
bool Conjugate,
int TriStorageOrder,
int OtherInnerStride,
83EIGEN_STRONG_INLINE
void trsmKernelR<Scalar,
Index, Mode, Conjugate, TriStorageOrder, OtherInnerStride,
84 Specialized>::kernel(
Index size,
Index otherSize,
const Scalar* _tri,
85 Index triStride, Scalar* _other,
Index otherIncr,
87 typedef typename NumTraits<Scalar>::Real RealScalar;
88 typedef blas_data_mapper<Scalar, Index, ColMajor, Unaligned, OtherInnerStride> LhsMapper;
89 typedef const_blas_data_mapper<Scalar, Index, TriStorageOrder> RhsMapper;
90 LhsMapper lhs(_other, otherStride, otherIncr);
91 RhsMapper rhs(_tri, triStride);
93 enum { RhsStorageOrder = TriStorageOrder, IsLower = (Mode &
Lower) ==
Lower };
94 conj_if<Conjugate>
conj;
96 for (
Index k = 0; k < size; ++k) {
97 Index j = IsLower ? size - k - 1 : k;
99 typename LhsMapper::LinearMapper r = lhs.getLinearMapper(0, j);
100 for (
Index k3 = 0; k3 < k; ++k3) {
101 Scalar b =
conj(rhs(IsLower ? j + 1 + k3 : k3, j));
102 typename LhsMapper::LinearMapper a = lhs.getLinearMapper(0, IsLower ? j + 1 + k3 : k3);
103 for (
Index i = 0; i < otherSize; ++i) r(i) -= a(i) * b;
106 Scalar inv_rjj = RealScalar(1) /
conj(rhs(j, j));
107 for (
Index i = 0; i < otherSize; ++i) r(i) *= inv_rjj;
113template <
typename Scalar,
typename Index,
int Side,
int Mode,
bool Conjugate,
int TriStorageOrder,
114 int OtherInnerStride>
115struct triangular_solve_matrix<Scalar,
Index, Side, Mode, Conjugate, TriStorageOrder,
RowMajor, OtherInnerStride> {
116 static void run(
Index size,
Index cols,
const Scalar* tri,
Index triStride, Scalar* _other,
Index otherIncr,
117 Index otherStride, level3_blocking<Scalar, Scalar>& blocking) {
118 triangular_solve_matrix<
121 OtherInnerStride>::run(size, cols, tri, triStride, _other, otherIncr, otherStride, blocking);
127template <
typename Scalar,
typename Index,
int Mode,
bool Conjugate,
int TriStorageOrder,
int OtherInnerStr
ide>
128struct triangular_solve_matrix<Scalar,
Index,
OnTheLeft, Mode, Conjugate, TriStorageOrder,
ColMajor, OtherInnerStride> {
129 static EIGEN_DONT_INLINE
void run(
Index size,
Index otherSize,
const Scalar* _tri,
Index triStride, Scalar* _other,
130 Index otherIncr,
Index otherStride, level3_blocking<Scalar, Scalar>& blocking);
133template <
typename Scalar,
typename Index,
int Mode,
bool Conjugate,
int TriStorageOrder,
int OtherInnerStr
ide>
134EIGEN_DONT_INLINE
void triangular_solve_matrix<Scalar,
Index,
OnTheLeft, Mode, Conjugate, TriStorageOrder,
ColMajor,
135 OtherInnerStride>::run(
Index size,
Index otherSize,
const Scalar* _tri,
136 Index triStride, Scalar* _other,
Index otherIncr,
138 level3_blocking<Scalar, Scalar>& blocking) {
139 Index cols = otherSize;
141 std::ptrdiff_t l1, l2, l3;
142 manage_caching_sizes(GetAction, &l1, &l2, &l3);
144#if defined(EIGEN_VECTORIZE_AVX512) && EIGEN_USE_AVX512_TRSM_L_KERNELS && EIGEN_ENABLE_AVX512_NOCOPY_TRSM_L_CUTOFFS
146 (OtherInnerStride == 1 && (std::is_same<Scalar, float>::value || std::is_same<Scalar, double>::value))) {
151 if (size < avx512_trsm_cutoff<Scalar>(l2, cols, L2Cap)) {
152 trsmKernelL<Scalar,
Index, Mode, Conjugate, TriStorageOrder, 1,
true>::kernel(
153 size, cols, _tri, triStride, _other, 1, otherStride);
159 typedef const_blas_data_mapper<Scalar, Index, TriStorageOrder> TriMapper;
160 typedef blas_data_mapper<Scalar, Index, ColMajor, Unaligned, OtherInnerStride> OtherMapper;
161 TriMapper tri(_tri, triStride);
162 OtherMapper other(_other, otherStride, otherIncr);
164 typedef gebp_traits<Scalar, Scalar> Traits;
166 enum { SmallPanelWidth = plain_enum_max(Traits::mr, Traits::nr), IsLower = (Mode &
Lower) ==
Lower };
168 Index kc = blocking.kc();
169 Index mc = (std::min)(size, blocking.mc());
171 std::size_t sizeA = kc * mc;
172 std::size_t sizeB = kc * cols;
174 ei_declare_aligned_stack_constructed_variable(Scalar, blockA, sizeA, blocking.blockA());
175 ei_declare_aligned_stack_constructed_variable(Scalar, blockB, sizeB, blocking.blockB());
177 gebp_kernel<Scalar, Scalar, Index, OtherMapper, Traits::mr, Traits::nr, Conjugate, false> gebp_kernel;
178 gemm_pack_lhs<Scalar,
Index, TriMapper, Traits::mr, Traits::LhsProgress,
typename Traits::LhsPacket4Packing,
181 gemm_pack_rhs<Scalar, Index, OtherMapper, Traits::nr, ColMajor, false, true> pack_rhs;
185 Index subcols = cols > 0 ? l2 / (4 *
sizeof(Scalar) * std::max<Index>(otherStride, size)) : 0;
186 subcols = std::max<Index>((subcols / Traits::nr) * Traits::nr, Traits::nr);
188 for (
Index k2 = IsLower ? 0 : size; IsLower ? k2 < size : k2 > 0; IsLower ? k2 += kc : k2 -= kc) {
189 const Index actual_kc = (std::min)(IsLower ? size - k2 : k2, kc);
204 for (
Index j2 = 0; j2 < cols; j2 += subcols) {
205 Index actual_cols = (std::min)(cols - j2, subcols);
207 for (
Index k1 = 0; k1 < actual_kc; k1 += SmallPanelWidth) {
208 Index actualPanelWidth = std::min<Index>(actual_kc - k1, SmallPanelWidth);
211 Index i = IsLower ? k2 + k1 : k2 - k1;
212#if defined(EIGEN_VECTORIZE_AVX512) && EIGEN_USE_AVX512_TRSM_L_KERNELS
214 (OtherInnerStride == 1 && (std::is_same<Scalar, float>::value || std::is_same<Scalar, double>::value))) {
215 i = IsLower ? k2 + k1 : k2 - k1 - actualPanelWidth;
218 trsmKernelL<Scalar,
Index, Mode, Conjugate, TriStorageOrder, OtherInnerStride,
true>::kernel(
219 actualPanelWidth, actual_cols, _tri + i + (i)*triStride, triStride,
220 _other + i * OtherInnerStride + j2 * otherStride, otherIncr, otherStride);
223 Index lengthTarget = actual_kc - k1 - actualPanelWidth;
224 Index startBlock = IsLower ? k2 + k1 : k2 - k1 - actualPanelWidth;
225 Index blockBOffset = IsLower ? k1 : lengthTarget;
228 pack_rhs(blockB + actual_kc * j2, other.getSubMapper(startBlock, j2), actualPanelWidth, actual_cols, actual_kc,
232 if (lengthTarget > 0) {
233 Index startTarget = IsLower ? k2 + k1 + actualPanelWidth : k2 - actual_kc;
235 pack_lhs(blockA, tri.getSubMapper(startTarget, startBlock), actualPanelWidth, lengthTarget);
237 gebp_kernel(other.getSubMapper(startTarget, j2), blockA, blockB + actual_kc * j2, lengthTarget,
238 actualPanelWidth, actual_cols, Scalar(-1), actualPanelWidth, actual_kc, 0, blockBOffset);
245 Index start = IsLower ? k2 + kc : 0;
246 Index end = IsLower ? size : k2 - kc;
247 for (
Index i2 = start; i2 < end; i2 += mc) {
248 const Index actual_mc = (std::min)(mc, end - i2);
250 pack_lhs(blockA, tri.getSubMapper(i2, IsLower ? k2 : k2 - kc), actual_kc, actual_mc);
252 gebp_kernel(other.getSubMapper(i2, 0), blockA, blockB, actual_mc, actual_kc, cols, Scalar(-1), -1, -1, 0, 0);
261template <
typename Scalar,
typename Index,
int Mode,
bool Conjugate,
int TriStorageOrder,
int OtherInnerStr
ide>
264 static EIGEN_DONT_INLINE
void run(
Index size,
Index otherSize,
const Scalar* _tri,
Index triStride, Scalar* _other,
265 Index otherIncr,
Index otherStride, level3_blocking<Scalar, Scalar>& blocking);
268template <
typename Scalar,
typename Index,
int Mode,
bool Conjugate,
int TriStorageOrder,
int OtherInnerStr
ide>
269EIGEN_DONT_INLINE
void triangular_solve_matrix<Scalar,
Index,
OnTheRight, Mode, Conjugate, TriStorageOrder,
ColMajor,
270 OtherInnerStride>::run(
Index size,
Index otherSize,
const Scalar* _tri,
271 Index triStride, Scalar* _other,
Index otherIncr,
273 level3_blocking<Scalar, Scalar>& blocking) {
274 Index rows = otherSize;
276#if defined(EIGEN_VECTORIZE_AVX512) && EIGEN_USE_AVX512_TRSM_R_KERNELS && EIGEN_ENABLE_AVX512_NOCOPY_TRSM_R_CUTOFFS
278 (OtherInnerStride == 1 && (std::is_same<Scalar, float>::value || std::is_same<Scalar, double>::value))) {
280 std::ptrdiff_t l1, l2, l3;
281 manage_caching_sizes(GetAction, &l1, &l2, &l3);
283 if (size < avx512_trsm_cutoff<Scalar>(l2, rows, L2Cap)) {
284 trsmKernelR<Scalar,
Index, Mode, Conjugate, TriStorageOrder, OtherInnerStride,
true>::kernel(
285 size, rows, _tri, triStride, _other, 1, otherStride);
291 typedef blas_data_mapper<Scalar, Index, ColMajor, Unaligned, OtherInnerStride> LhsMapper;
292 typedef const_blas_data_mapper<Scalar, Index, TriStorageOrder> RhsMapper;
293 LhsMapper lhs(_other, otherStride, otherIncr);
294 RhsMapper rhs(_tri, triStride);
296 typedef gebp_traits<Scalar, Scalar> Traits;
298 RhsStorageOrder = TriStorageOrder,
299 SmallPanelWidth = plain_enum_max(Traits::mr, Traits::nr),
303 Index kc = blocking.kc();
304 Index mc = (std::min)(rows, blocking.mc());
306 std::size_t sizeA = kc * mc;
307 std::size_t sizeB = kc * size;
309 ei_declare_aligned_stack_constructed_variable(Scalar, blockA, sizeA, blocking.blockA());
310 ei_declare_aligned_stack_constructed_variable(Scalar, blockB, sizeB, blocking.blockB());
312 gebp_kernel<Scalar, Scalar, Index, LhsMapper, Traits::mr, Traits::nr, false, Conjugate> gebp_kernel;
313 gemm_pack_rhs<Scalar, Index, RhsMapper, Traits::nr, RhsStorageOrder> pack_rhs;
314 gemm_pack_rhs<Scalar, Index, RhsMapper, Traits::nr, RhsStorageOrder, false, true> pack_rhs_panel;
315 gemm_pack_lhs<Scalar,
Index, LhsMapper, Traits::mr, Traits::LhsProgress,
typename Traits::LhsPacket4Packing,
ColMajor,
319 for (
Index k2 = IsLower ? size : 0; IsLower ? k2 > 0 : k2 < size; IsLower ? k2 -= kc : k2 += kc) {
320 const Index actual_kc = (std::min)(IsLower ? k2 : size - k2, kc);
321 Index actual_k2 = IsLower ? k2 - actual_kc : k2;
323 Index startPanel = IsLower ? 0 : k2 + actual_kc;
324 Index rs = IsLower ? actual_k2 : size - actual_k2 - actual_kc;
325 Scalar* geb = blockB + actual_kc * actual_kc;
327 if (rs > 0) pack_rhs(geb, rhs.getSubMapper(actual_k2, startPanel), actual_kc, rs);
332 for (
Index j2 = 0; j2 < actual_kc; j2 += SmallPanelWidth) {
333 Index actualPanelWidth = std::min<Index>(actual_kc - j2, SmallPanelWidth);
334 Index actual_j2 = actual_k2 + j2;
335 Index panelOffset = IsLower ? j2 + actualPanelWidth : 0;
336 Index panelLength = IsLower ? actual_kc - j2 - actualPanelWidth : j2;
339 pack_rhs_panel(blockB + j2 * actual_kc, rhs.getSubMapper(actual_k2 + panelOffset, actual_j2), panelLength,
340 actualPanelWidth, actual_kc, panelOffset);
344 for (
Index i2 = 0; i2 < rows; i2 += mc) {
345 const Index actual_mc = (std::min)(mc, rows - i2);
350 for (
Index j2 = IsLower ? (actual_kc - ((actual_kc % SmallPanelWidth) ?
Index(actual_kc % SmallPanelWidth)
351 :
Index(SmallPanelWidth)))
353 IsLower ? j2 >= 0 : j2 < actual_kc; IsLower ? j2 -= SmallPanelWidth : j2 += SmallPanelWidth) {
354 Index actualPanelWidth = std::min<Index>(actual_kc - j2, SmallPanelWidth);
355 Index absolute_j2 = actual_k2 + j2;
356 Index panelOffset = IsLower ? j2 + actualPanelWidth : 0;
357 Index panelLength = IsLower ? actual_kc - j2 - actualPanelWidth : j2;
360 if (panelLength > 0) {
361 gebp_kernel(lhs.getSubMapper(i2, absolute_j2), blockA, blockB + j2 * actual_kc, actual_mc, panelLength,
362 actualPanelWidth, Scalar(-1), actual_kc, actual_kc,
363 panelOffset, panelOffset);
368 trsmKernelR<Scalar,
Index, Mode, Conjugate, TriStorageOrder, OtherInnerStride,
369 true>::kernel(actualPanelWidth, actual_mc,
370 _tri + absolute_j2 + absolute_j2 * triStride, triStride,
371 _other + i2 * OtherInnerStride + absolute_j2 * otherStride,
372 otherIncr, otherStride);
375 pack_lhs_panel(blockA, lhs.getSubMapper(i2, absolute_j2), actualPanelWidth, actual_mc, actual_kc, j2);
380 gebp_kernel(lhs.getSubMapper(i2, startPanel), blockA, geb, actual_mc, actual_kc, rs, Scalar(-1), -1, -1, 0, 0);
@ UnitDiag
Definition Constants.h:215
@ Lower
Definition Constants.h:211
@ Upper
Definition Constants.h:213
@ ColMajor
Definition Constants.h:318
@ RowMajor
Definition Constants.h:320
@ OnTheLeft
Definition Constants.h:331
@ OnTheRight
Definition Constants.h:333
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82