level2_cplx_impl.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009-2010 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 #include "common.h"
11 
19 EIGEN_BLAS_FUNC(hemv)
20 (const char *uplo, const int *n, const RealScalar *palpha, const RealScalar *pa, const int *lda, const RealScalar *px,
21  const int *incx, const RealScalar *pbeta, RealScalar *py, const int *incy) {
22  typedef void (*functype)(int, const Scalar *, int, const Scalar *, Scalar *, Scalar);
23  static const functype func[2] = {
24  // array index: UP
26  false>::run),
27  // array index: LO
29  false>::run),
30  };
31 
32  const Scalar *a = reinterpret_cast<const Scalar *>(pa);
33  const Scalar *x = reinterpret_cast<const Scalar *>(px);
34  Scalar *y = reinterpret_cast<Scalar *>(py);
35  Scalar alpha = *reinterpret_cast<const Scalar *>(palpha);
36  Scalar beta = *reinterpret_cast<const Scalar *>(pbeta);
37 
38  // check arguments
39  int info = 0;
40  if (UPLO(*uplo) == INVALID)
41  info = 1;
42  else if (*n < 0)
43  info = 2;
44  else if (*lda < std::max(1, *n))
45  info = 5;
46  else if (*incx == 0)
47  info = 7;
48  else if (*incy == 0)
49  info = 10;
50  if (info) return xerbla_(SCALAR_SUFFIX_UP "HEMV ", &info);
51 
52  if (*n == 0) return;
53 
56 
57  if (beta != Scalar(1)) {
58  if (beta == Scalar(0))
59  make_vector(actual_y, *n).setZero();
60  else
61  make_vector(actual_y, *n) *= beta;
62  }
63 
64  if (alpha != Scalar(0)) {
65  int code = UPLO(*uplo);
66  if (code >= 2 || func[code] == 0) return;
67 
68  func[code](*n, a, *lda, actual_x, actual_y, alpha);
69  }
70 
71  if (actual_x != x) delete[] actual_x;
72  if (actual_y != y) delete[] copy_back(actual_y, y, *n, *incy);
73 }
74 
82 // EIGEN_BLAS_FUNC(hbmv)(char *uplo, int *n, int *k, RealScalar *alpha, RealScalar *a, int *lda,
83 // RealScalar *x, int *incx, RealScalar *beta, RealScalar *y, int *incy)
84 // {
85 // return 1;
86 // }
87 
95 // EIGEN_BLAS_FUNC(hpmv)(char *uplo, int *n, RealScalar *alpha, RealScalar *ap, RealScalar *x, int *incx, RealScalar
96 // *beta, RealScalar *y, int *incy)
97 // {
98 // return 1;
99 // }
100 
109  typedef void (*functype)(int, Scalar *, const Scalar *, RealScalar);
110  static const functype func[2] = {
111  // array index: UP
113  // array index: LO
115  };
116 
117  Scalar *x = reinterpret_cast<Scalar *>(px);
118  Scalar *ap = reinterpret_cast<Scalar *>(pap);
120 
121  int info = 0;
122  if (UPLO(*uplo) == INVALID)
123  info = 1;
124  else if (*n < 0)
125  info = 2;
126  else if (*incx == 0)
127  info = 5;
128  if (info) return xerbla_(SCALAR_SUFFIX_UP "HPR ", &info);
129 
130  if (alpha == Scalar(0)) return;
131 
133 
134  int code = UPLO(*uplo);
135  if (code >= 2 || func[code] == 0) return;
136 
137  func[code](*n, ap, x_cpy, alpha);
138 
139  if (x_cpy != x) delete[] x_cpy;
140 }
141 
149 EIGEN_BLAS_FUNC(hpr2)
150 (char *uplo, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pap) {
151  typedef void (*functype)(int, Scalar *, const Scalar *, const Scalar *, Scalar);
152  static const functype func[2] = {
153  // array index: UP
155  // array index: LO
157  };
158 
159  Scalar *x = reinterpret_cast<Scalar *>(px);
160  Scalar *y = reinterpret_cast<Scalar *>(py);
161  Scalar *ap = reinterpret_cast<Scalar *>(pap);
162  Scalar alpha = *reinterpret_cast<Scalar *>(palpha);
163 
164  int info = 0;
165  if (UPLO(*uplo) == INVALID)
166  info = 1;
167  else if (*n < 0)
168  info = 2;
169  else if (*incx == 0)
170  info = 5;
171  else if (*incy == 0)
172  info = 7;
173  if (info) return xerbla_(SCALAR_SUFFIX_UP "HPR2 ", &info);
174 
175  if (alpha == Scalar(0)) return;
176 
179 
180  int code = UPLO(*uplo);
181  if (code >= 2 || func[code] == 0) return;
182 
183  func[code](*n, ap, x_cpy, y_cpy, alpha);
184 
185  if (x_cpy != x) delete[] x_cpy;
186  if (y_cpy != y) delete[] y_cpy;
187 }
188 
197  typedef void (*functype)(int, Scalar *, int, const Scalar *, const Scalar *, const Scalar &);
198  static const functype func[2] = {
199  // array index: UP
201  // array index: LO
203  };
204 
205  Scalar *x = reinterpret_cast<Scalar *>(px);
206  Scalar *a = reinterpret_cast<Scalar *>(pa);
207  RealScalar alpha = *reinterpret_cast<RealScalar *>(palpha);
208 
209  int info = 0;
210  if (UPLO(*uplo) == INVALID)
211  info = 1;
212  else if (*n < 0)
213  info = 2;
214  else if (*incx == 0)
215  info = 5;
216  else if (*lda < std::max(1, *n))
217  info = 7;
218  if (info) return xerbla_(SCALAR_SUFFIX_UP "HER ", &info);
219 
220  if (alpha == RealScalar(0)) return;
221 
223 
224  int code = UPLO(*uplo);
225  if (code >= 2 || func[code] == 0) return;
226 
227  func[code](*n, a, *lda, x_cpy, x_cpy, alpha);
228 
229  matrix(a, *n, *n, *lda).diagonal().imag().setZero();
230 
231  if (x_cpy != x) delete[] x_cpy;
232 }
233 
241 EIGEN_BLAS_FUNC(her2)
242 (char *uplo, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pa,
243  int *lda) {
244  typedef void (*functype)(int, Scalar *, int, const Scalar *, const Scalar *, Scalar);
245  static const functype func[2] = {
246  // array index: UP
248  // array index: LO
250  };
251 
252  Scalar *x = reinterpret_cast<Scalar *>(px);
253  Scalar *y = reinterpret_cast<Scalar *>(py);
254  Scalar *a = reinterpret_cast<Scalar *>(pa);
255  Scalar alpha = *reinterpret_cast<Scalar *>(palpha);
256 
257  int info = 0;
258  if (UPLO(*uplo) == INVALID)
259  info = 1;
260  else if (*n < 0)
261  info = 2;
262  else if (*incx == 0)
263  info = 5;
264  else if (*incy == 0)
265  info = 7;
266  else if (*lda < std::max(1, *n))
267  info = 9;
268  if (info) return xerbla_(SCALAR_SUFFIX_UP "HER2 ", &info);
269 
270  if (alpha == Scalar(0)) return;
271 
274 
275  int code = UPLO(*uplo);
276  if (code >= 2 || func[code] == 0) return;
277 
278  func[code](*n, a, *lda, x_cpy, y_cpy, alpha);
279 
280  matrix(a, *n, *n, *lda).diagonal().imag().setZero();
281 
282  if (x_cpy != x) delete[] x_cpy;
283  if (y_cpy != y) delete[] y_cpy;
284 }
285 
293 EIGEN_BLAS_FUNC(geru)
294 (int *m, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pa, int *lda) {
295  Scalar *x = reinterpret_cast<Scalar *>(px);
296  Scalar *y = reinterpret_cast<Scalar *>(py);
297  Scalar *a = reinterpret_cast<Scalar *>(pa);
298  Scalar alpha = *reinterpret_cast<Scalar *>(palpha);
299 
300  int info = 0;
301  if (*m < 0)
302  info = 1;
303  else if (*n < 0)
304  info = 2;
305  else if (*incx == 0)
306  info = 5;
307  else if (*incy == 0)
308  info = 7;
309  else if (*lda < std::max(1, *m))
310  info = 9;
311  if (info) return xerbla_(SCALAR_SUFFIX_UP "GERU ", &info);
312 
313  if (alpha == Scalar(0)) return;
314 
317 
319  alpha);
320 
321  if (x_cpy != x) delete[] x_cpy;
322  if (y_cpy != y) delete[] y_cpy;
323 }
324 
332 EIGEN_BLAS_FUNC(gerc)
333 (int *m, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pa, int *lda) {
334  Scalar *x = reinterpret_cast<Scalar *>(px);
335  Scalar *y = reinterpret_cast<Scalar *>(py);
336  Scalar *a = reinterpret_cast<Scalar *>(pa);
337  Scalar alpha = *reinterpret_cast<Scalar *>(palpha);
338 
339  int info = 0;
340  if (*m < 0)
341  info = 1;
342  else if (*n < 0)
343  info = 2;
344  else if (*incx == 0)
345  info = 5;
346  else if (*incy == 0)
347  info = 7;
348  else if (*lda < std::max(1, *m))
349  info = 9;
350  if (info) return xerbla_(SCALAR_SUFFIX_UP "GERC ", &info);
351 
352  if (alpha == Scalar(0)) return;
353 
356 
358  alpha);
359 
360  if (x_cpy != x) delete[] x_cpy;
361  if (y_cpy != y) delete[] y_cpy;
362 }
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
#define SCALAR_SUFFIX_UP
Definition: blas/complex_double.cpp:12
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, 1 >, 0, Eigen::InnerStride< Eigen::Dynamic > > make_vector(T *data, int size, int incr)
Definition: common.h:99
T * copy_back(T *x_cpy, T *x, int n, int incx)
Definition: common.h:136
T * get_compact_vector(T *x, int n, int incx)
Definition: common.h:124
#define INVALID
Definition: common.h:52
#define UPLO(X)
Definition: common.h:59
#define max(a, b)
Definition: datatypes.h:23
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ ColMajor
Definition: Constants.h:318
return int(ret)+1
EIGEN_BLAS_FUNC(EIGEN_CAT(REAL_SCALAR_SUFFIX, scal))(int *n
EIGEN_BLAS_FUNC() hpr(char *uplo, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *pap)
Definition: level2_cplx_impl.h:108
const char const int const RealScalar const RealScalar const int const RealScalar const int const RealScalar RealScalar * py
Definition: level2_cplx_impl.h:21
const char * uplo
Definition: level2_cplx_impl.h:20
const char const int const RealScalar const RealScalar const int const RealScalar * px
Definition: level2_cplx_impl.h:20
const Scalar * x
Definition: level2_cplx_impl.h:33
Scalar * x_cpy
Definition: level2_cplx_impl.h:177
Scalar * y_cpy
Definition: level2_cplx_impl.h:178
const Scalar * a
Definition: level2_cplx_impl.h:32
const char const int * n
Definition: level2_cplx_impl.h:20
int code
Definition: level2_cplx_impl.h:180
const char const int const RealScalar const RealScalar const int const RealScalar const int const RealScalar * pbeta
Definition: level2_cplx_impl.h:21
const char const int const RealScalar const RealScalar const int * lda
Definition: level2_cplx_impl.h:20
const char const int const RealScalar const RealScalar * pa
Definition: level2_cplx_impl.h:20
Scalar * y
Definition: level2_cplx_impl.h:34
int * m
Definition: level2_cplx_impl.h:294
int info
Definition: level2_cplx_impl.h:39
Scalar * ap
Definition: level2_cplx_impl.h:161
EIGEN_BLAS_FUNC() her(char *uplo, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *pa, int *lda)
Definition: level2_cplx_impl.h:196
Scalar * actual_y
Definition: level2_cplx_impl.h:55
Scalar alpha
Definition: level2_cplx_impl.h:35
const char const int const RealScalar * palpha
Definition: level2_cplx_impl.h:20
char int RealScalar RealScalar int RealScalar int RealScalar * pap
Definition: level2_cplx_impl.h:150
matrix(a, *n, *n, *lda).diagonal().imag().setZero()
Scalar beta
Definition: level2_cplx_impl.h:36
const Scalar * actual_x
Definition: level2_cplx_impl.h:54
const char const int const RealScalar const RealScalar const int const RealScalar const int const RealScalar RealScalar const int * incy
Definition: level2_cplx_impl.h:21
const char const int const RealScalar const RealScalar const int const RealScalar const int * incx
Definition: level2_cplx_impl.h:21
void(* functype)(DenseIndex, DenseIndex, DenseIndex, const Scalar *, DenseIndex, const Scalar *, DenseIndex, Scalar *, DenseIndex, DenseIndex, Scalar, Eigen::internal::level3_blocking< Scalar, Scalar > &, Eigen::internal::GemmParallelInfo< DenseIndex > *)
Definition: level3_impl.h:22
Definition: Rank2Update.h:20
Definition: SelfadjointMatrixVector.h:34
Definition: PackedSelfadjointProduct.h:20
Definition: GeneralMatrixMatrixTriangular.h:19
Definition: benchGeometry.cpp:21
void run(const string &dir_name, LinearSolver *linear_solver_pt, const unsigned nel_1d, bool mess_up_order)
Definition: two_d_poisson_compare_solvers.cc:317
EIGEN_WEAK_LINKING void xerbla_(const char *msg, int *info)
Definition: xerbla.cpp:14