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