// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2009-2010 Gael Guennebaud // // This Source Code Form is subject to the terms of the Mozilla // Public License v. 2.0. If a copy of the MPL was not distributed // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. #include "common.h" template struct general_matrix_vector_product_wrapper { static void run(Index rows, Index cols,const Scalar *lhs, Index lhsStride, const Scalar *rhs, Index rhsIncr, Scalar* res, Index resIncr, Scalar alpha) { typedef internal::const_blas_data_mapper LhsMapper; typedef internal::const_blas_data_mapper RhsMapper; internal::general_matrix_vector_product ::run( rows, cols, LhsMapper(lhs, lhsStride), RhsMapper(rhs, rhsIncr), res, resIncr, alpha); } }; int EIGEN_BLAS_FUNC(gemv)(const char *opa, const int *m, const int *n, const RealScalar *palpha, const RealScalar *pa, const int *lda, const RealScalar *pb, const int *incb, const RealScalar *pbeta, RealScalar *pc, const int *incc) { typedef void (*functype)(int, int, const Scalar *, int, const Scalar *, int , Scalar *, int, Scalar); static const functype func[4] = { // array index: NOTR (general_matrix_vector_product_wrapper::run), // array index: TR (general_matrix_vector_product_wrapper::run), // array index: ADJ (general_matrix_vector_product_wrapper::run), 0 }; const Scalar* a = reinterpret_cast(pa); const Scalar* b = reinterpret_cast(pb); Scalar* c = reinterpret_cast(pc); Scalar alpha = *reinterpret_cast(palpha); Scalar beta = *reinterpret_cast(pbeta); // check arguments int info = 0; if(OP(*opa)==INVALID) info = 1; else if(*m<0) info = 2; else if(*n<0) info = 3; else if(*lda=4 || func[code]==0) return 0; func[code](actual_m, actual_n, a, *lda, actual_b, 1, actual_c, 1, alpha); if(actual_b!=b) delete[] actual_b; if(actual_c!=c) delete[] copy_back(actual_c,c,actual_m,*incc); return 1; } int EIGEN_BLAS_FUNC(trsv)(const char *uplo, const char *opa, const char *diag, const int *n, const RealScalar *pa, const int *lda, RealScalar *pb, const int *incb) { typedef void (*functype)(int, const Scalar *, int, Scalar *); static const functype func[16] = { // array index: NOTR | (UP << 2) | (NUNIT << 3) (internal::triangular_solve_vector::run), // array index: TR | (UP << 2) | (NUNIT << 3) (internal::triangular_solve_vector::run), // array index: ADJ | (UP << 2) | (NUNIT << 3) (internal::triangular_solve_vector::run), 0, // array index: NOTR | (LO << 2) | (NUNIT << 3) (internal::triangular_solve_vector::run), // array index: TR | (LO << 2) | (NUNIT << 3) (internal::triangular_solve_vector::run), // array index: ADJ | (LO << 2) | (NUNIT << 3) (internal::triangular_solve_vector::run), 0, // array index: NOTR | (UP << 2) | (UNIT << 3) (internal::triangular_solve_vector::run), // array index: TR | (UP << 2) | (UNIT << 3) (internal::triangular_solve_vector::run), // array index: ADJ | (UP << 2) | (UNIT << 3) (internal::triangular_solve_vector::run), 0, // array index: NOTR | (LO << 2) | (UNIT << 3) (internal::triangular_solve_vector::run), // array index: TR | (LO << 2) | (UNIT << 3) (internal::triangular_solve_vector::run), // array index: ADJ | (LO << 2) | (UNIT << 3) (internal::triangular_solve_vector::run), 0 }; const Scalar* a = reinterpret_cast(pa); Scalar* b = reinterpret_cast(pb); int info = 0; if(UPLO(*uplo)==INVALID) info = 1; else if(OP(*opa)==INVALID) info = 2; else if(DIAG(*diag)==INVALID) info = 3; else if(*n<0) info = 4; else if(*lda::run), // array index: TR | (UP << 2) | (NUNIT << 3) (internal::triangular_matrix_vector_product::run), // array index: ADJ | (UP << 2) | (NUNIT << 3) (internal::triangular_matrix_vector_product::run), 0, // array index: NOTR | (LO << 2) | (NUNIT << 3) (internal::triangular_matrix_vector_product::run), // array index: TR | (LO << 2) | (NUNIT << 3) (internal::triangular_matrix_vector_product::run), // array index: ADJ | (LO << 2) | (NUNIT << 3) (internal::triangular_matrix_vector_product::run), 0, // array index: NOTR | (UP << 2) | (UNIT << 3) (internal::triangular_matrix_vector_product::run), // array index: TR | (UP << 2) | (UNIT << 3) (internal::triangular_matrix_vector_product::run), // array index: ADJ | (UP << 2) | (UNIT << 3) (internal::triangular_matrix_vector_product::run), 0, // array index: NOTR | (LO << 2) | (UNIT << 3) (internal::triangular_matrix_vector_product::run), // array index: TR | (LO << 2) | (UNIT << 3) (internal::triangular_matrix_vector_product::run), // array index: ADJ | (LO << 2) | (UNIT << 3) (internal::triangular_matrix_vector_product::run), 0 }; const Scalar* a = reinterpret_cast(pa); Scalar* b = reinterpret_cast(pb); int info = 0; if(UPLO(*uplo)==INVALID) info = 1; else if(OP(*opa)==INVALID) info = 2; else if(DIAG(*diag)==INVALID) info = 3; else if(*n<0) info = 4; else if(*lda res(*n); res.setZero(); int code = OP(*opa) | (UPLO(*uplo) << 2) | (DIAG(*diag) << 3); if(code>=16 || func[code]==0) return 0; func[code](*n, *n, a, *lda, actual_b, 1, res.data(), 1, Scalar(1)); copy_back(res.data(),b,*n,*incb); if(actual_b!=b) delete[] actual_b; return 1; } /** GBMV performs one of the matrix-vector operations * * y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y, * * where alpha and beta are scalars, x and y are vectors and A is an * m by n band matrix, with kl sub-diagonals and ku super-diagonals. */ int EIGEN_BLAS_FUNC(gbmv)(char *trans, int *m, int *n, int *kl, int *ku, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *px, int *incx, RealScalar *pbeta, RealScalar *py, int *incy) { const Scalar* a = reinterpret_cast(pa); const Scalar* x = reinterpret_cast(px); Scalar* y = reinterpret_cast(py); Scalar alpha = *reinterpret_cast(palpha); Scalar beta = *reinterpret_cast(pbeta); int coeff_rows = *kl+*ku+1; int info = 0; if(OP(*trans)==INVALID) info = 1; else if(*m<0) info = 2; else if(*n<0) info = 3; else if(*kl<0) info = 4; else if(*ku<0) info = 5; else if(*lda(pa); Scalar* x = reinterpret_cast(px); int coeff_rows = *k + 1; int info = 0; if(UPLO(*uplo)==INVALID) info = 1; else if(OP(*opa)==INVALID) info = 2; else if(DIAG(*diag)==INVALID) info = 3; else if(*n<0) info = 4; else if(*k<0) info = 5; else if(*lda::run), // array index: TR | (UP << 2) | (NUNIT << 3) (internal::band_solve_triangular_selector::run), // array index: ADJ | (UP << 2) | (NUNIT << 3) (internal::band_solve_triangular_selector::run), 0, // array index: NOTR | (LO << 2) | (NUNIT << 3) (internal::band_solve_triangular_selector::run), // array index: TR | (LO << 2) | (NUNIT << 3) (internal::band_solve_triangular_selector::run), // array index: ADJ | (LO << 2) | (NUNIT << 3) (internal::band_solve_triangular_selector::run), 0, // array index: NOTR | (UP << 2) | (UNIT << 3) (internal::band_solve_triangular_selector::run), // array index: TR | (UP << 2) | (UNIT << 3) (internal::band_solve_triangular_selector::run), // array index: ADJ | (UP << 2) | (UNIT << 3) (internal::band_solve_triangular_selector::run), 0, // array index: NOTR | (LO << 2) | (UNIT << 3) (internal::band_solve_triangular_selector::run), // array index: TR | (LO << 2) | (UNIT << 3) (internal::band_solve_triangular_selector::run), // array index: ADJ | (LO << 2) | (UNIT << 3) (internal::band_solve_triangular_selector::run), 0, }; Scalar* a = reinterpret_cast(pa); Scalar* x = reinterpret_cast(px); int coeff_rows = *k+1; int info = 0; if(UPLO(*uplo)==INVALID) info = 1; else if(OP(*op)==INVALID) info = 2; else if(DIAG(*diag)==INVALID) info = 3; else if(*n<0) info = 4; else if(*k<0) info = 5; else if(*lda=16 || func[code]==0) return 0; func[code](*n, *k, a, *lda, actual_x); if(actual_x!=x) delete[] copy_back(actual_x,x,actual_n,*incx); return 0; } /** DTPMV performs one of the matrix-vector operations * * x := A*x, or x := A'*x, * * where x is an n element vector and A is an n by n unit, or non-unit, * upper or lower triangular matrix, supplied in packed form. */ int EIGEN_BLAS_FUNC(tpmv)(char *uplo, char *opa, char *diag, int *n, RealScalar *pap, RealScalar *px, int *incx) { typedef void (*functype)(int, const Scalar*, const Scalar*, Scalar*, Scalar); static const functype func[16] = { // array index: NOTR | (UP << 2) | (NUNIT << 3) (internal::packed_triangular_matrix_vector_product::run), // array index: TR | (UP << 2) | (NUNIT << 3) (internal::packed_triangular_matrix_vector_product::run), // array index: ADJ | (UP << 2) | (NUNIT << 3) (internal::packed_triangular_matrix_vector_product::run), 0, // array index: NOTR | (LO << 2) | (NUNIT << 3) (internal::packed_triangular_matrix_vector_product::run), // array index: TR | (LO << 2) | (NUNIT << 3) (internal::packed_triangular_matrix_vector_product::run), // array index: ADJ | (LO << 2) | (NUNIT << 3) (internal::packed_triangular_matrix_vector_product::run), 0, // array index: NOTR | (UP << 2) | (UNIT << 3) (internal::packed_triangular_matrix_vector_product::run), // array index: TR | (UP << 2) | (UNIT << 3) (internal::packed_triangular_matrix_vector_product::run), // array index: ADJ | (UP << 2) | (UNIT << 3) (internal::packed_triangular_matrix_vector_product::run), 0, // array index: NOTR | (LO << 2) | (UNIT << 3) (internal::packed_triangular_matrix_vector_product::run), // array index: TR | (LO << 2) | (UNIT << 3) (internal::packed_triangular_matrix_vector_product::run), // array index: ADJ | (LO << 2) | (UNIT << 3) (internal::packed_triangular_matrix_vector_product::run), 0 }; Scalar* ap = reinterpret_cast(pap); Scalar* x = reinterpret_cast(px); int info = 0; if(UPLO(*uplo)==INVALID) info = 1; else if(OP(*opa)==INVALID) info = 2; else if(DIAG(*diag)==INVALID) info = 3; else if(*n<0) info = 4; else if(*incx==0) info = 7; if(info) return xerbla_(SCALAR_SUFFIX_UP"TPMV ",&info,6); if(*n==0) return 1; Scalar* actual_x = get_compact_vector(x,*n,*incx); Matrix res(*n); res.setZero(); int code = OP(*opa) | (UPLO(*uplo) << 2) | (DIAG(*diag) << 3); if(code>=16 || func[code]==0) return 0; func[code](*n, ap, actual_x, res.data(), Scalar(1)); copy_back(res.data(),x,*n,*incx); if(actual_x!=x) delete[] actual_x; return 1; } /** DTPSV solves one of the systems of equations * * A*x = b, or A'*x = b, * * where b and x are n element vectors and A is an n by n unit, or * non-unit, upper or lower triangular matrix, supplied in packed form. * * No test for singularity or near-singularity is included in this * routine. Such tests must be performed before calling this routine. */ int EIGEN_BLAS_FUNC(tpsv)(char *uplo, char *opa, char *diag, int *n, RealScalar *pap, RealScalar *px, int *incx) { typedef void (*functype)(int, const Scalar*, Scalar*); static const functype func[16] = { // array index: NOTR | (UP << 2) | (NUNIT << 3) (internal::packed_triangular_solve_vector::run), // array index: TR | (UP << 2) | (NUNIT << 3) (internal::packed_triangular_solve_vector::run), // array index: ADJ | (UP << 2) | (NUNIT << 3) (internal::packed_triangular_solve_vector::run), 0, // array index: NOTR | (LO << 2) | (NUNIT << 3) (internal::packed_triangular_solve_vector::run), // array index: TR | (LO << 2) | (NUNIT << 3) (internal::packed_triangular_solve_vector::run), // array index: ADJ | (LO << 2) | (NUNIT << 3) (internal::packed_triangular_solve_vector::run), 0, // array index: NOTR | (UP << 2) | (UNIT << 3) (internal::packed_triangular_solve_vector::run), // array index: TR | (UP << 2) | (UNIT << 3) (internal::packed_triangular_solve_vector::run), // array index: ADJ | (UP << 2) | (UNIT << 3) (internal::packed_triangular_solve_vector::run), 0, // array index: NOTR | (LO << 2) | (UNIT << 3) (internal::packed_triangular_solve_vector::run), // array index: TR | (LO << 2) | (UNIT << 3) (internal::packed_triangular_solve_vector::run), // array index: ADJ | (LO << 2) | (UNIT << 3) (internal::packed_triangular_solve_vector::run), 0 }; Scalar* ap = reinterpret_cast(pap); Scalar* x = reinterpret_cast(px); int info = 0; if(UPLO(*uplo)==INVALID) info = 1; else if(OP(*opa)==INVALID) info = 2; else if(DIAG(*diag)==INVALID) info = 3; else if(*n<0) info = 4; else if(*incx==0) info = 7; if(info) return xerbla_(SCALAR_SUFFIX_UP"TPSV ",&info,6); Scalar* actual_x = get_compact_vector(x,*n,*incx); int code = OP(*opa) | (UPLO(*uplo) << 2) | (DIAG(*diag) << 3); func[code](*n, ap, actual_x); if(actual_x!=x) delete[] copy_back(actual_x,x,*n,*incx); return 1; }