Eigen  5.0.1-dev+60122df6
 
Loading...
Searching...
No Matches
SelfadjointRank2Update.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_SELFADJOINTRANK2UPTADE_H
11#define EIGEN_SELFADJOINTRANK2UPTADE_H
12
13// IWYU pragma: private
14#include "../InternalHeaderCheck.h"
15
16namespace Eigen {
17
18namespace internal {
19
20/* Optimized selfadjoint matrix += alpha * uv' + conj(alpha)*vu'
21 * It corresponds to the Level2 syr2 BLAS routine
22 */
23
24template <typename Scalar, typename Index, typename UType, typename VType, int UpLo>
25struct selfadjoint_rank2_update_selector;
26
27template <typename Scalar, typename Index, typename UType, typename VType>
28struct selfadjoint_rank2_update_selector<Scalar, Index, UType, VType, Lower> {
29 static EIGEN_DEVICE_FUNC void run(Scalar* mat, Index stride, const UType& u, const VType& v, const Scalar& alpha) {
30 const Index size = u.size();
31 for (Index i = 0; i < size; ++i) {
32 Map<Matrix<Scalar, Dynamic, 1>>(mat + stride * i + i, size - i) +=
33 (numext::conj(alpha) * numext::conj(u.coeff(i))) * v.tail(size - i) +
34 (alpha * numext::conj(v.coeff(i))) * u.tail(size - i);
35 }
36 }
37};
38
39template <typename Scalar, typename Index, typename UType, typename VType>
40struct selfadjoint_rank2_update_selector<Scalar, Index, UType, VType, Upper> {
41 static void run(Scalar* mat, Index stride, const UType& u, const VType& v, const Scalar& alpha) {
42 const Index size = u.size();
43 for (Index i = 0; i < size; ++i)
44 Map<Matrix<Scalar, Dynamic, 1>>(mat + stride * i, i + 1) +=
45 (numext::conj(alpha) * numext::conj(u.coeff(i))) * v.head(i + 1) +
46 (alpha * numext::conj(v.coeff(i))) * u.head(i + 1);
47 }
48};
49
50template <bool Cond, typename T>
51using conj_expr_if =
52 std::conditional<!Cond, const T&, CwiseUnaryOp<scalar_conjugate_op<typename traits<T>::Scalar>, T>>;
53
54} // end namespace internal
55
56template <typename MatrixType, unsigned int UpLo>
57template <typename DerivedU, typename DerivedV>
59 const MatrixBase<DerivedU>& u, const MatrixBase<DerivedV>& v, const Scalar& alpha) {
60 typedef internal::blas_traits<DerivedU> UBlasTraits;
61 typedef typename UBlasTraits::DirectLinearAccessType ActualUType;
62 typedef internal::remove_all_t<ActualUType> ActualUType_;
63 internal::add_const_on_value_type_t<ActualUType> actualU = UBlasTraits::extract(u.derived());
64
65 typedef internal::blas_traits<DerivedV> VBlasTraits;
66 typedef typename VBlasTraits::DirectLinearAccessType ActualVType;
67 typedef internal::remove_all_t<ActualVType> ActualVType_;
68 internal::add_const_on_value_type_t<ActualVType> actualV = VBlasTraits::extract(v.derived());
69
70 // If MatrixType is row major, then we use the routine for lower triangular in the upper triangular case and
71 // vice versa, and take the complex conjugate of all coefficients and vector entries.
72
73 enum { IsRowMajor = (internal::traits<MatrixType>::Flags & RowMajorBit) ? 1 : 0 };
74 Scalar actualAlpha = alpha * UBlasTraits::extractScalarFactor(u.derived()) *
75 numext::conj(VBlasTraits::extractScalarFactor(v.derived()));
76 if (IsRowMajor) actualAlpha = numext::conj(actualAlpha);
77
78 typedef internal::remove_all_t<
79 typename internal::conj_expr_if<int(IsRowMajor) ^ int(UBlasTraits::NeedToConjugate), ActualUType_>::type>
80 UType;
81 typedef internal::remove_all_t<
82 typename internal::conj_expr_if<int(IsRowMajor) ^ int(VBlasTraits::NeedToConjugate), ActualVType_>::type>
83 VType;
84 internal::selfadjoint_rank2_update_selector<Scalar, Index, UType, VType,
85 (IsRowMajor ? int(UpLo == Upper ? Lower : Upper)
86 : UpLo)>::run(_expression().const_cast_derived().data(),
87 _expression().outerStride(), UType(actualU),
88 VType(actualV), actualAlpha);
89
90 return *this;
91}
92
93} // end namespace Eigen
94
95#endif // EIGEN_SELFADJOINTRANK2UPTADE_H
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
Expression of a selfadjoint matrix from a triangular part of a dense matrix.
Definition SelfAdjointView.h:51
SelfAdjointView & rankUpdate(const MatrixBase< DerivedU > &u, const MatrixBase< DerivedV > &v, const Scalar &alpha=Scalar(1))
@ Lower
Definition Constants.h:211
@ Upper
Definition Constants.h:213
const unsigned int RowMajorBit
Definition Constants.h:70
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82