level2_real_impl.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // Eigen is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 3 of the License, or (at your option) any later version.
00010 //
00011 // Alternatively, you can redistribute it and/or
00012 // modify it under the terms of the GNU General Public License as
00013 // published by the Free Software Foundation; either version 2 of
00014 // the License, or (at your option) any later version.
00015 //
00016 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00017 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00018 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00019 // GNU General Public License for more details.
00020 //
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License and a copy of the GNU General Public License along with
00023 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00024 
00025 #include "common.h"
00026 
00027 // y = alpha*A*x + beta*y
00028 int EIGEN_BLAS_FUNC(symv) (char *uplo, int *n, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *px, int *incx, RealScalar *pbeta, RealScalar *py, int *incy)
00029 {
00030   Scalar* a = reinterpret_cast<Scalar*>(pa);
00031   Scalar* x = reinterpret_cast<Scalar*>(px);
00032   Scalar* y = reinterpret_cast<Scalar*>(py);
00033   Scalar alpha  = *reinterpret_cast<Scalar*>(palpha);
00034   Scalar beta   = *reinterpret_cast<Scalar*>(pbeta);
00035 
00036   // check arguments
00037   int info = 0;
00038   if(UPLO(*uplo)==INVALID)        info = 1;
00039   else if(*n<0)                   info = 2;
00040   else if(*lda<std::max(1,*n))    info = 5;
00041   else if(*incx==0)               info = 7;
00042   else if(*incy==0)               info = 10;
00043   if(info)
00044     return xerbla_(SCALAR_SUFFIX_UP"SYMV ",&info,6);
00045 
00046   if(*n==0)
00047     return 0;
00048 
00049   Scalar* actual_x = get_compact_vector(x,*n,*incx);
00050   Scalar* actual_y = get_compact_vector(y,*n,*incy);
00051 
00052   if(beta!=Scalar(1))
00053   {
00054     if(beta==Scalar(0)) vector(actual_y, *n).setZero();
00055     else                vector(actual_y, *n) *= beta;
00056   }
00057 
00058   // TODO performs a direct call to the underlying implementation function
00059        if(UPLO(*uplo)==UP) vector(actual_y,*n).noalias() += matrix(a,*n,*n,*lda).selfadjointView<Upper>() * (alpha * vector(actual_x,*n));
00060   else if(UPLO(*uplo)==LO) vector(actual_y,*n).noalias() += matrix(a,*n,*n,*lda).selfadjointView<Lower>() * (alpha * vector(actual_x,*n));
00061 
00062   if(actual_x!=x) delete[] actual_x;
00063   if(actual_y!=y) delete[] copy_back(actual_y,y,*n,*incy);
00064 
00065   return 1;
00066 }
00067 
00068 // C := alpha*x*x' + C
00069 int EIGEN_BLAS_FUNC(syr)(char *uplo, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *pc, int *ldc)
00070 {
00071 
00072 //   typedef void (*functype)(int, const Scalar *, int, Scalar *, int, Scalar);
00073 //   static functype func[2];
00074 
00075 //   static bool init = false;
00076 //   if(!init)
00077 //   {
00078 //     for(int k=0; k<2; ++k)
00079 //       func[k] = 0;
00080 //
00081 //     func[UP] = (internal::selfadjoint_product<Scalar,ColMajor,ColMajor,false,UpperTriangular>::run);
00082 //     func[LO] = (internal::selfadjoint_product<Scalar,ColMajor,ColMajor,false,LowerTriangular>::run);
00083 
00084 //     init = true;
00085 //   }
00086 
00087   Scalar* x = reinterpret_cast<Scalar*>(px);
00088   Scalar* c = reinterpret_cast<Scalar*>(pc);
00089   Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
00090 
00091   int info = 0;
00092   if(UPLO(*uplo)==INVALID)                                            info = 1;
00093   else if(*n<0)                                                       info = 2;
00094   else if(*incx==0)                                                   info = 5;
00095   else if(*ldc<std::max(1,*n))                                        info = 7;
00096   if(info)
00097     return xerbla_(SCALAR_SUFFIX_UP"SYR  ",&info,6);
00098 
00099   if(*n==0 || alpha==Scalar(0)) return 1;
00100 
00101   // if the increment is not 1, let's copy it to a temporary vector to enable vectorization
00102   Scalar* x_cpy = get_compact_vector(x,*n,*incx);
00103 
00104   Matrix<Scalar,Dynamic,Dynamic> m2(matrix(c,*n,*n,*ldc));
00105   
00106   // TODO check why this is not accurate enough for lapack tests
00107 //   if(UPLO(*uplo)==LO)       matrix(c,*n,*n,*ldc).selfadjointView<Lower>().rankUpdate(vector(x_cpy,*n), alpha);
00108 //   else if(UPLO(*uplo)==UP)  matrix(c,*n,*n,*ldc).selfadjointView<Upper>().rankUpdate(vector(x_cpy,*n), alpha);
00109 
00110   if(UPLO(*uplo)==LO)
00111     for(int j=0;j<*n;++j)
00112       matrix(c,*n,*n,*ldc).col(j).tail(*n-j) += alpha * x_cpy[j] * vector(x_cpy+j,*n-j);
00113   else
00114     for(int j=0;j<*n;++j)
00115       matrix(c,*n,*n,*ldc).col(j).head(j+1) += alpha * x_cpy[j] * vector(x_cpy,j+1);
00116 
00117   if(x_cpy!=x)  delete[] x_cpy;
00118 
00119   return 1;
00120 }
00121 
00122 // C := alpha*x*y' + alpha*y*x' + C
00123 int EIGEN_BLAS_FUNC(syr2)(char *uplo, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pc, int *ldc)
00124 {
00125 //   typedef void (*functype)(int, const Scalar *, int, const Scalar *, int, Scalar *, int, Scalar);
00126 //   static functype func[2];
00127 //
00128 //   static bool init = false;
00129 //   if(!init)
00130 //   {
00131 //     for(int k=0; k<2; ++k)
00132 //       func[k] = 0;
00133 //
00134 //     func[UP] = (internal::selfadjoint_product<Scalar,ColMajor,ColMajor,false,UpperTriangular>::run);
00135 //     func[LO] = (internal::selfadjoint_product<Scalar,ColMajor,ColMajor,false,LowerTriangular>::run);
00136 //
00137 //     init = true;
00138 //   }
00139 
00140   Scalar* x = reinterpret_cast<Scalar*>(px);
00141   Scalar* y = reinterpret_cast<Scalar*>(py);
00142   Scalar* c = reinterpret_cast<Scalar*>(pc);
00143   Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
00144 
00145   int info = 0;
00146   if(UPLO(*uplo)==INVALID)                                            info = 1;
00147   else if(*n<0)                                                       info = 2;
00148   else if(*incx==0)                                                   info = 5;
00149   else if(*incy==0)                                                   info = 7;
00150   else if(*ldc<std::max(1,*n))                                        info = 9;
00151   if(info)
00152     return xerbla_(SCALAR_SUFFIX_UP"SYR2 ",&info,6);
00153 
00154   if(alpha==Scalar(0))
00155     return 1;
00156 
00157   Scalar* x_cpy = get_compact_vector(x,*n,*incx);
00158   Scalar* y_cpy = get_compact_vector(y,*n,*incy);
00159 
00160   // TODO perform direct calls to underlying implementation
00161   if(UPLO(*uplo)==LO)       matrix(c,*n,*n,*ldc).selfadjointView<Lower>().rankUpdate(vector(x_cpy,*n), vector(y_cpy,*n), alpha);
00162   else if(UPLO(*uplo)==UP)  matrix(c,*n,*n,*ldc).selfadjointView<Upper>().rankUpdate(vector(x_cpy,*n), vector(y_cpy,*n), alpha);
00163 
00164   if(x_cpy!=x)  delete[] x_cpy;
00165   if(y_cpy!=y)  delete[] y_cpy;
00166 
00167 //   int code = UPLO(*uplo);
00168 //   if(code>=2 || func[code]==0)
00169 //     return 0;
00170 
00171 //   func[code](*n, a, *inca, b, *incb, c, *ldc, alpha);
00172   return 1;
00173 }
00174 
00182 // int EIGEN_BLAS_FUNC(sbmv)( char *uplo, int *n, int *k, RealScalar *alpha, RealScalar *a, int *lda,
00183 //                            RealScalar *x, int *incx, RealScalar *beta, RealScalar *y, int *incy)
00184 // {
00185 //   return 1;
00186 // }
00187 
00188 
00197 // int EIGEN_BLAS_FUNC(spmv)(char *uplo, int *n, RealScalar *alpha, RealScalar *ap, RealScalar *x, int *incx, RealScalar *beta, RealScalar *y, int *incy)
00198 // {
00199 //   return 1;
00200 // }
00201 
00209 // int EIGEN_BLAS_FUNC(spr)(char *uplo, int *n, Scalar *alpha, Scalar *x, int *incx, Scalar *ap)
00210 // {
00211 //   return 1;
00212 // }
00213 
00221 // int EIGEN_BLAS_FUNC(spr2)(char *uplo, int *n, RealScalar *alpha, RealScalar *x, int *incx, RealScalar *y, int *incy, RealScalar *ap)
00222 // {
00223 //   return 1;
00224 // }
00225 


re_vision
Author(s): Dorian Galvez-Lopez
autogenerated on Sun Jan 5 2014 11:31:38