<T>LAPACK 0.1.1
C++ Template Linear Algebra PACKage
Loading...
Searching...
No Matches
her.hpp
Go to the documentation of this file.
1
3//
4// Copyright (c) 2017-2021, University of Tennessee. All rights reserved.
5// Copyright (c) 2021-2023, University of Colorado Denver. All rights reserved.
6//
7// This file is part of <T>LAPACK.
8// <T>LAPACK is free software: you can redistribute it and/or modify it under
9// the terms of the BSD 3-Clause license. See the accompanying LICENSE file.
10
11#ifndef TLAPACK_BLAS_HER_HH
12#define TLAPACK_BLAS_HER_HH
13
15
16namespace tlapack {
17
42template <TLAPACK_MATRIX matrixA_t,
43 TLAPACK_VECTOR vectorX_t,
44 TLAPACK_REAL alpha_t,
45 enable_if_t<(
46 /* Requires: */
47 is_real<alpha_t>),
48 int> = 0,
49 class T = type_t<matrixA_t>,
50 disable_if_allow_optblas_t<pair<alpha_t, real_type<T> >,
51 pair<matrixA_t, T>,
52 pair<vectorX_t, T> > = 0>
53void her(Uplo uplo, const alpha_t& alpha, const vectorX_t& x, matrixA_t& A)
54{
55 // data traits
56 using idx_t = size_type<matrixA_t>;
58
59 // constants
60 const idx_t n = nrows(A);
61
62 // check arguments
63 tlapack_check_false(uplo != Uplo::Lower && uplo != Uplo::Upper);
64 tlapack_check_false(size(x) != n);
65 tlapack_check_false(ncols(A) != n);
66
67 if (uplo == Uplo::Upper) {
68 for (idx_t j = 0; j < n; ++j) {
69 const scalar_t tmp = alpha * conj(x[j]);
70 for (idx_t i = 0; i < j; ++i)
71 A(i, j) += x[i] * tmp;
72 A(j, j) =
73 real(A(j, j)) + real(x[j]) * real(tmp) - imag(x[j]) * imag(tmp);
74 }
75 }
76 else {
77 for (idx_t j = 0; j < n; ++j) {
78 const scalar_t tmp = alpha * conj(x[j]);
79 A(j, j) =
80 real(A(j, j)) + real(x[j]) * real(tmp) - imag(x[j]) * imag(tmp);
81 for (idx_t i = j + 1; i < n; ++i)
82 A(i, j) += x[i] * tmp;
83 }
84 }
85}
86
87#ifdef TLAPACK_USE_LAPACKPP
88
89template <TLAPACK_LEGACY_MATRIX matrixA_t,
90 TLAPACK_LEGACY_VECTOR vectorX_t,
91 TLAPACK_REAL alpha_t,
92 class T = type_t<matrixA_t>,
93 enable_if_allow_optblas_t<pair<alpha_t, real_type<T> >,
94 pair<matrixA_t, T>,
95 pair<vectorX_t, T> > = 0>
96void her(Uplo uplo, const alpha_t alpha, const vectorX_t& x, matrixA_t& A)
97{
98 // Legacy objects
99 auto A_ = legacy_matrix(A);
100 auto x_ = legacy_vector(x);
101
102 // Constants to forward
103 constexpr Layout L = layout<matrixA_t>;
104 const auto& n = A_.n;
105
106 return ::blas::her((::blas::Layout)L, (::blas::Uplo)uplo, n, alpha, x_.ptr,
107 x_.inc, A_.ptr, A_.ldim);
108}
109
110#endif
111
112} // namespace tlapack
113
114#endif // #ifndef TLAPACK_BLAS_HER_HH
Uplo
Definition types.hpp:45
constexpr real_type< T > real(const T &x) noexcept
Extends std::real() to real datatypes.
Definition utils.hpp:71
constexpr T conj(const T &x) noexcept
Extends std::conj() to real datatypes.
Definition utils.hpp:100
constexpr real_type< T > imag(const T &x) noexcept
Extends std::imag() to real datatypes.
Definition utils.hpp:86
#define TLAPACK_LEGACY_VECTOR
Macro for tlapack::concepts::LegacyVector compatible with C++17.
Definition concepts.hpp:954
#define TLAPACK_LEGACY_MATRIX
Macro for tlapack::concepts::LegacyMatrix compatible with C++17.
Definition concepts.hpp:951
#define TLAPACK_VECTOR
Macro for tlapack::concepts::Vector compatible with C++17.
Definition concepts.hpp:906
#define TLAPACK_MATRIX
Macro for tlapack::concepts::Matrix compatible with C++17.
Definition concepts.hpp:896
#define TLAPACK_REAL
Macro for tlapack::concepts::Real compatible with C++17.
Definition concepts.hpp:918
void her(Uplo uplo, const alpha_t &alpha, const vectorX_t &x, matrixA_t &A)
Hermitian matrix rank-1 update:
Definition her.hpp:53
#define tlapack_check_false(cond)
Throw an error if cond is true.
Definition exceptionHandling.hpp:113
Concept for types that represent tlapack::Uplo.
typename traits::real_type_traits< Types..., int >::type real_type
The common real type of the list of types.
Definition scalar_type_traits.hpp:113