1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2009-2015 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_BLAS_COMMON_H
11#define EIGEN_BLAS_COMMON_H
12
13#include "../Eigen/Core"
14#include "../Eigen/Jacobi"
15
16#include <complex>
17
18#ifndef SCALAR
19#error the token SCALAR must be defined to compile this file
20#endif
21
22#include "../Eigen/src/misc/blas.h"
23
24#define NOTR    0
25#define TR      1
26#define ADJ     2
27
28#define LEFT    0
29#define RIGHT   1
30
31#define UP      0
32#define LO      1
33
34#define NUNIT   0
35#define UNIT    1
36
37#define INVALID 0xff
38
39#define OP(X)   (   ((X)=='N' || (X)=='n') ? NOTR   \
40                  : ((X)=='T' || (X)=='t') ? TR     \
41                  : ((X)=='C' || (X)=='c') ? ADJ    \
42                  : INVALID)
43
44#define SIDE(X) (   ((X)=='L' || (X)=='l') ? LEFT   \
45                  : ((X)=='R' || (X)=='r') ? RIGHT  \
46                  : INVALID)
47
48#define UPLO(X) (   ((X)=='U' || (X)=='u') ? UP     \
49                  : ((X)=='L' || (X)=='l') ? LO     \
50                  : INVALID)
51
52#define DIAG(X) (   ((X)=='N' || (X)=='n') ? NUNIT  \
53                  : ((X)=='U' || (X)=='u') ? UNIT   \
54                  : INVALID)
55
56
57inline bool check_op(const char* op)
58{
59  return OP(*op)!=0xff;
60}
61
62inline bool check_side(const char* side)
63{
64  return SIDE(*side)!=0xff;
65}
66
67inline bool check_uplo(const char* uplo)
68{
69  return UPLO(*uplo)!=0xff;
70}
71
72
73namespace Eigen {
74#include "BandTriangularSolver.h"
75#include "GeneralRank1Update.h"
76#include "PackedSelfadjointProduct.h"
77#include "PackedTriangularMatrixVector.h"
78#include "PackedTriangularSolverVector.h"
79#include "Rank2Update.h"
80}
81
82using namespace Eigen;
83
84typedef SCALAR Scalar;
85typedef NumTraits<Scalar>::Real RealScalar;
86typedef std::complex<RealScalar> Complex;
87
88enum
89{
90  IsComplex = Eigen::NumTraits<SCALAR>::IsComplex,
91  Conj = IsComplex
92};
93
94typedef Matrix<Scalar,Dynamic,Dynamic,ColMajor> PlainMatrixType;
95typedef Map<Matrix<Scalar,Dynamic,Dynamic,ColMajor>, 0, OuterStride<> > MatrixType;
96typedef Map<const Matrix<Scalar,Dynamic,Dynamic,ColMajor>, 0, OuterStride<> > ConstMatrixType;
97typedef Map<Matrix<Scalar,Dynamic,1>, 0, InnerStride<Dynamic> > StridedVectorType;
98typedef Map<Matrix<Scalar,Dynamic,1> > CompactVectorType;
99
100template<typename T>
101Map<Matrix<T,Dynamic,Dynamic,ColMajor>, 0, OuterStride<> >
102matrix(T* data, int rows, int cols, int stride)
103{
104  return Map<Matrix<T,Dynamic,Dynamic,ColMajor>, 0, OuterStride<> >(data, rows, cols, OuterStride<>(stride));
105}
106
107template<typename T>
108Map<const Matrix<T,Dynamic,Dynamic,ColMajor>, 0, OuterStride<> >
109matrix(const T* data, int rows, int cols, int stride)
110{
111  return Map<const Matrix<T,Dynamic,Dynamic,ColMajor>, 0, OuterStride<> >(data, rows, cols, OuterStride<>(stride));
112}
113
114template<typename T>
115Map<Matrix<T,Dynamic,1>, 0, InnerStride<Dynamic> > make_vector(T* data, int size, int incr)
116{
117  return Map<Matrix<T,Dynamic,1>, 0, InnerStride<Dynamic> >(data, size, InnerStride<Dynamic>(incr));
118}
119
120template<typename T>
121Map<const Matrix<T,Dynamic,1>, 0, InnerStride<Dynamic> > make_vector(const T* data, int size, int incr)
122{
123  return Map<const Matrix<T,Dynamic,1>, 0, InnerStride<Dynamic> >(data, size, InnerStride<Dynamic>(incr));
124}
125
126template<typename T>
127Map<Matrix<T,Dynamic,1> > make_vector(T* data, int size)
128{
129  return Map<Matrix<T,Dynamic,1> >(data, size);
130}
131
132template<typename T>
133Map<const Matrix<T,Dynamic,1> > make_vector(const T* data, int size)
134{
135  return Map<const Matrix<T,Dynamic,1> >(data, size);
136}
137
138template<typename T>
139T* get_compact_vector(T* x, int n, int incx)
140{
141  if(incx==1)
142    return x;
143
144  typename Eigen::internal::remove_const<T>::type* ret = new Scalar[n];
145  if(incx<0) make_vector(ret,n) = make_vector(x,n,-incx).reverse();
146  else       make_vector(ret,n) = make_vector(x,n, incx);
147  return ret;
148}
149
150template<typename T>
151T* copy_back(T* x_cpy, T* x, int n, int incx)
152{
153  if(x_cpy==x)
154    return 0;
155
156  if(incx<0) make_vector(x,n,-incx).reverse() = make_vector(x_cpy,n);
157  else       make_vector(x,n, incx)           = make_vector(x_cpy,n);
158  return x_cpy;
159}
160
161#define EIGEN_BLAS_FUNC(X) EIGEN_CAT(SCALAR_SUFFIX,X##_)
162
163#endif // EIGEN_BLAS_COMMON_H
164