<T>LAPACK 0.1.1
C++ Template Linear Algebra PACKage
Loading...
Searching...
No Matches
hemv.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_HEMV_HH
12#define TLAPACK_BLAS_HEMV_HH
13
15
16namespace tlapack {
17
42template <TLAPACK_MATRIX matrixA_t,
43 TLAPACK_VECTOR vectorX_t,
44 TLAPACK_VECTOR vectorY_t,
45 TLAPACK_SCALAR alpha_t,
46 TLAPACK_SCALAR beta_t,
47 class T = type_t<vectorY_t>,
48 disable_if_allow_optblas_t<pair<alpha_t, T>,
49 pair<matrixA_t, T>,
50 pair<vectorX_t, T>,
51 pair<vectorY_t, T>,
52 pair<beta_t, T> > = 0>
54 const alpha_t& alpha,
55 const matrixA_t& A,
56 const vectorX_t& x,
57 const beta_t& beta,
58 vectorY_t& y)
59{
60 // data traits
61 using TA = type_t<matrixA_t>;
62 using TX = type_t<vectorX_t>;
63 using idx_t = size_type<matrixA_t>;
64
65 // constants
66 const idx_t n = nrows(A);
67
68 // check arguments
69 tlapack_check_false(uplo != Uplo::Lower && uplo != Uplo::Upper);
70 tlapack_check_false(ncols(A) != n);
71 tlapack_check_false(size(x) != n);
72 tlapack_check_false(size(y) != n);
73
74 // form y = beta*y
75 for (idx_t i = 0; i < n; ++i)
76 y[i] *= beta;
77
78 if (uplo == Uplo::Upper) {
79 // A is stored in upper triangle
80 // form y += alpha * A * x
81 for (idx_t j = 0; j < n; ++j) {
84 for (idx_t i = 0; i < j; ++i) {
85 y[i] += tmp1 * A(i, j);
86 sum += conj(A(i, j)) * x[i];
87 }
88 y[j] += tmp1 * real(A(j, j)) + alpha * sum;
89 }
90 }
91 else {
92 // A is stored in lower triangle
93 // form y += alpha * A * x
94 for (idx_t j = 0; j < n; ++j) {
97 for (idx_t i = j + 1; i < n; ++i) {
98 y[i] += tmp1 * A(i, j);
99 sum += conj(A(i, j)) * x[i];
100 }
101 y[j] += tmp1 * real(A(j, j)) + alpha * sum;
102 }
103 }
104}
105
106#ifdef TLAPACK_USE_LAPACKPP
107
108template <TLAPACK_LEGACY_MATRIX matrixA_t,
109 TLAPACK_LEGACY_VECTOR vectorX_t,
110 TLAPACK_LEGACY_VECTOR vectorY_t,
111 TLAPACK_SCALAR alpha_t,
112 TLAPACK_SCALAR beta_t,
113 class T = type_t<vectorY_t>,
114 enable_if_allow_optblas_t<pair<alpha_t, T>,
115 pair<matrixA_t, T>,
116 pair<vectorX_t, T>,
117 pair<vectorY_t, T>,
118 pair<beta_t, T> > = 0>
119void hemv(Uplo uplo,
120 const alpha_t alpha,
121 const matrixA_t& A,
122 const vectorX_t& x,
123 const beta_t beta,
124 vectorY_t& y)
125{
126 // Legacy objects
127 auto A_ = legacy_matrix(A);
128 auto x_ = legacy_vector(x);
129 auto y_ = legacy_vector(y);
130
131 // Constants to forward
132 constexpr Layout L = layout<matrixA_t>;
133 const auto& n = A_.n;
134
135 // Warnings for NaNs and Infs
136 if (alpha == alpha_t(0))
138 -2, "Infs and NaNs in A or x will not propagate to y on output");
139 if (beta == beta_t(0) && !is_same_v<beta_t, StrongZero>)
141 -5,
142 "Infs and NaNs in y on input will not propagate to y on output");
143
144 return ::blas::hemv((::blas::Layout)L, (::blas::Uplo)uplo, n, alpha, A_.ptr,
145 A_.ldim, x_.ptr, x_.inc, (T)beta, y_.ptr, y_.inc);
146}
147
148#endif
149
173template <TLAPACK_MATRIX matrixA_t,
174 TLAPACK_VECTOR vectorX_t,
175 TLAPACK_VECTOR vectorY_t,
176 TLAPACK_SCALAR alpha_t>
178 const alpha_t& alpha,
179 const matrixA_t& A,
180 const vectorX_t& x,
181 vectorY_t& y)
182{
183 return hemv(uplo, alpha, A, x, StrongZero(), y);
184}
185
186} // namespace tlapack
187
188#endif // #ifndef TLAPACK_BLAS_HEMV_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
#define TLAPACK_SCALAR
Macro for tlapack::concepts::Scalar compatible with C++17.
Definition concepts.hpp:915
#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
void hemv(Uplo uplo, const alpha_t &alpha, const matrixA_t &A, const vectorX_t &x, const beta_t &beta, vectorY_t &y)
Hermitian matrix-vector multiply:
Definition hemv.hpp:53
#define tlapack_check_false(cond)
Throw an error if cond is true.
Definition exceptionHandling.hpp:113
#define tlapack_warning(info, detailedInfo)
Warning handler.
Definition exceptionHandling.hpp:156
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
Strong zero type.
Definition StrongZero.hpp:43