level2_real_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 
12 // y = alpha*A*x + beta*y
13 EIGEN_BLAS_FUNC(symv)
14 (const char *uplo, const int *n, const RealScalar *palpha, const RealScalar *pa, const int *lda, const RealScalar *px,
15  const int *incx, const RealScalar *pbeta, RealScalar *py, const int *incy) {
16  typedef void (*functype)(int, const Scalar *, int, const Scalar *, Scalar *, Scalar);
17  using Eigen::ColMajor;
18  using Eigen::Lower;
19  using Eigen::Upper;
20  static const functype func[2] = {
21  // array index: UP
23  // array index: LO
25  };
26 
27  const Scalar *a = reinterpret_cast<const Scalar *>(pa);
28  const Scalar *x = reinterpret_cast<const Scalar *>(px);
29  Scalar *y = reinterpret_cast<Scalar *>(py);
30  Scalar alpha = *reinterpret_cast<const Scalar *>(palpha);
31  Scalar beta = *reinterpret_cast<const Scalar *>(pbeta);
32 
33  // check arguments
34  int info = 0;
35  if (UPLO(*uplo) == INVALID)
36  info = 1;
37  else if (*n < 0)
38  info = 2;
39  else if (*lda < std::max(1, *n))
40  info = 5;
41  else if (*incx == 0)
42  info = 7;
43  else if (*incy == 0)
44  info = 10;
45  if (info) return xerbla_(SCALAR_SUFFIX_UP "SYMV ", &info);
46 
47  if (*n == 0) return;
48 
51 
52  if (beta != Scalar(1)) {
53  if (beta == Scalar(0))
54  make_vector(actual_y, *n).setZero();
55  else
56  make_vector(actual_y, *n) *= beta;
57  }
58 
59  int code = UPLO(*uplo);
60  if (code >= 2 || func[code] == 0) return;
61 
62  func[code](*n, a, *lda, actual_x, actual_y, alpha);
63 
64  if (actual_x != x) delete[] actual_x;
65  if (actual_y != y) delete[] copy_back(actual_y, y, *n, *incy);
66 }
67 
68 // C := alpha*x*x' + C
69 EIGEN_BLAS_FUNC(syr)
70 (const char *uplo, const int *n, const RealScalar *palpha, const RealScalar *px, const int *incx, RealScalar *pc,
71  const int *ldc) {
72  typedef void (*functype)(int, Scalar *, int, const Scalar *, const Scalar *, const Scalar &);
73  using Eigen::ColMajor;
74  using Eigen::Lower;
75  using Eigen::Upper;
76  static const functype func[2] = {
77  // array index: UP
79  // array index: LO
81  };
82 
83  const Scalar *x = reinterpret_cast<const Scalar *>(px);
84  Scalar *c = reinterpret_cast<Scalar *>(pc);
85  Scalar alpha = *reinterpret_cast<const Scalar *>(palpha);
86 
87  int info = 0;
88  if (UPLO(*uplo) == INVALID)
89  info = 1;
90  else if (*n < 0)
91  info = 2;
92  else if (*incx == 0)
93  info = 5;
94  else if (*ldc < std::max(1, *n))
95  info = 7;
96  if (info) return xerbla_(SCALAR_SUFFIX_UP "SYR ", &info);
97 
98  if (*n == 0 || alpha == Scalar(0)) return;
99 
100  // if the increment is not 1, let's copy it to a temporary vector to enable vectorization
101  const Scalar *x_cpy = get_compact_vector(x, *n, *incx);
102 
103  int code = UPLO(*uplo);
104  if (code >= 2 || func[code] == 0) return;
105 
106  func[code](*n, c, *ldc, x_cpy, x_cpy, alpha);
107 
108  if (x_cpy != x) delete[] x_cpy;
109 }
110 
111 // C := alpha*x*y' + alpha*y*x' + C
112 EIGEN_BLAS_FUNC(syr2)
113 (const char *uplo, const int *n, const RealScalar *palpha, const RealScalar *px, const int *incx, const RealScalar *py,
114  const int *incy, RealScalar *pc, const int *ldc) {
115  typedef void (*functype)(int, Scalar *, int, const Scalar *, const Scalar *, Scalar);
116  static const functype func[2] = {
117  // array index: UP
119  // array index: LO
121  };
122 
123  const Scalar *x = reinterpret_cast<const Scalar *>(px);
124  const Scalar *y = reinterpret_cast<const Scalar *>(py);
125  Scalar *c = reinterpret_cast<Scalar *>(pc);
126  Scalar alpha = *reinterpret_cast<const Scalar *>(palpha);
127 
128  int info = 0;
129  if (UPLO(*uplo) == INVALID)
130  info = 1;
131  else if (*n < 0)
132  info = 2;
133  else if (*incx == 0)
134  info = 5;
135  else if (*incy == 0)
136  info = 7;
137  else if (*ldc < std::max(1, *n))
138  info = 9;
139  if (info) return xerbla_(SCALAR_SUFFIX_UP "SYR2 ", &info);
140 
141  if (alpha == Scalar(0)) return;
142 
145 
146  int code = UPLO(*uplo);
147  if (code >= 2 || func[code] == 0) return;
148 
149  func[code](*n, c, *ldc, x_cpy, y_cpy, alpha);
150 
151  if (x_cpy != x) delete[] x_cpy;
152  if (y_cpy != y) delete[] y_cpy;
153 
154  // int code = UPLO(*uplo);
155  // if(code>=2 || func[code]==0)
156  // return 0;
157 
158  // func[code](*n, a, *inca, b, *incb, c, *ldc, alpha);
159 }
160 
168 // EIGEN_BLAS_FUNC(sbmv)( char *uplo, int *n, int *k, RealScalar *alpha, RealScalar *a, int *lda,
169 // RealScalar *x, int *incx, RealScalar *beta, RealScalar *y, int *incy)
170 // {
171 // return 1;
172 // }
173 
182 // EIGEN_BLAS_FUNC(spmv)(char *uplo, int *n, RealScalar *alpha, RealScalar *ap, RealScalar *x, int *incx, RealScalar
183 // *beta, RealScalar *y, int *incy)
184 // {
185 // return 1;
186 // }
187 
195 EIGEN_BLAS_FUNC(spr)(char *uplo, int *n, Scalar *palpha, Scalar *px, int *incx, Scalar *pap) {
196  typedef void (*functype)(int, Scalar *, const Scalar *, Scalar);
197  static const functype func[2] = {
198  // array index: UP
200  // array index: LO
202  };
203 
204  Scalar *x = reinterpret_cast<Scalar *>(px);
205  Scalar *ap = reinterpret_cast<Scalar *>(pap);
206  Scalar alpha = *reinterpret_cast<Scalar *>(palpha);
207 
208  int info = 0;
209  if (UPLO(*uplo) == INVALID)
210  info = 1;
211  else if (*n < 0)
212  info = 2;
213  else if (*incx == 0)
214  info = 5;
215  if (info) return xerbla_(SCALAR_SUFFIX_UP "SPR ", &info);
216 
217  if (alpha == Scalar(0)) return;
218 
220 
221  int code = UPLO(*uplo);
222  if (code >= 2 || func[code] == 0) return;
223 
224  func[code](*n, ap, x_cpy, alpha);
225 
226  if (x_cpy != x) delete[] x_cpy;
227 }
228 
236 EIGEN_BLAS_FUNC(spr2)
237 (char *uplo, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pap) {
238  typedef void (*functype)(int, Scalar *, const Scalar *, const Scalar *, Scalar);
239  static const functype func[2] = {
240  // array index: UP
242  // array index: LO
244  };
245 
246  Scalar *x = reinterpret_cast<Scalar *>(px);
247  Scalar *y = reinterpret_cast<Scalar *>(py);
248  Scalar *ap = reinterpret_cast<Scalar *>(pap);
249  Scalar alpha = *reinterpret_cast<Scalar *>(palpha);
250 
251  int info = 0;
252  if (UPLO(*uplo) == INVALID)
253  info = 1;
254  else if (*n < 0)
255  info = 2;
256  else if (*incx == 0)
257  info = 5;
258  else if (*incy == 0)
259  info = 7;
260  if (info) return xerbla_(SCALAR_SUFFIX_UP "SPR2 ", &info);
261 
262  if (alpha == Scalar(0)) return;
263 
266 
267  int code = UPLO(*uplo);
268  if (code >= 2 || func[code] == 0) return;
269 
270  func[code](*n, ap, x_cpy, y_cpy, alpha);
271 
272  if (x_cpy != x) delete[] x_cpy;
273  if (y_cpy != y) delete[] y_cpy;
274 }
275 
283 EIGEN_BLAS_FUNC(ger)
284 (int *m, int *n, Scalar *palpha, Scalar *px, int *incx, Scalar *py, int *incy, Scalar *pa, int *lda) {
285  Scalar *x = reinterpret_cast<Scalar *>(px);
286  Scalar *y = reinterpret_cast<Scalar *>(py);
287  Scalar *a = reinterpret_cast<Scalar *>(pa);
288  Scalar alpha = *reinterpret_cast<Scalar *>(palpha);
289 
290  int info = 0;
291  if (*m < 0)
292  info = 1;
293  else if (*n < 0)
294  info = 2;
295  else if (*incx == 0)
296  info = 5;
297  else if (*incy == 0)
298  info = 7;
299  else if (*lda < std::max(1, *m))
300  info = 9;
301  if (info) return xerbla_(SCALAR_SUFFIX_UP "GER ", &info);
302 
303  if (alpha == Scalar(0)) return;
304 
307 
309  alpha);
310 
311  if (x_cpy != x) delete[] x_cpy;
312  if (y_cpy != y) delete[] y_cpy;
313 }
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
const char const int const RealScalar const RealScalar const int const RealScalar const int const RealScalar RealScalar const int * incy
Definition: level2_real_impl.h:15
const char const int const RealScalar * palpha
Definition: level2_real_impl.h:14
const char const int const RealScalar const RealScalar const int * lda
Definition: level2_real_impl.h:14
const Scalar * y_cpy
Definition: level2_real_impl.h:144
const Scalar * x_cpy
Definition: level2_real_impl.h:143
const Scalar * x
Definition: level2_real_impl.h:28
const char const int const RealScalar const RealScalar const int const RealScalar const int const RealScalar * pbeta
Definition: level2_real_impl.h:15
const char const int const RealScalar const RealScalar const int const RealScalar * px
Definition: level2_real_impl.h:14
const char const int * n
Definition: level2_real_impl.h:14
int code
Definition: level2_real_impl.h:59
const Scalar * a
Definition: level2_real_impl.h:27
Scalar * y
Definition: level2_real_impl.h:29
int info
Definition: level2_real_impl.h:34
const char * uplo
Definition: level2_real_impl.h:14
EIGEN_BLAS_FUNC() spr(char *uplo, int *n, Scalar *palpha, Scalar *px, int *incx, Scalar *pap)
Definition: level2_real_impl.h:195
Scalar * ap
Definition: level2_real_impl.h:248
Scalar * actual_y
Definition: level2_real_impl.h:50
const char const int const RealScalar const RealScalar const int const RealScalar const int * incx
Definition: level2_real_impl.h:15
Scalar alpha
Definition: level2_real_impl.h:30
Scalar beta
Definition: level2_real_impl.h:31
Scalar * c
Definition: level2_real_impl.h:84
char int RealScalar RealScalar int RealScalar int RealScalar * pap
Definition: level2_real_impl.h:237
const char const int const RealScalar const RealScalar const int RealScalar const int * ldc
Definition: level2_real_impl.h:71
const char const int const RealScalar const RealScalar const int RealScalar * pc
Definition: level2_real_impl.h:70
const Scalar * actual_x
Definition: level2_real_impl.h:49
const char const int const RealScalar const RealScalar * pa
Definition: level2_real_impl.h:14
const char const int const RealScalar const RealScalar const int const RealScalar const int const RealScalar RealScalar * py
Definition: level2_real_impl.h:15
int * m
Definition: level2_real_impl.h:284
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