// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr> // // 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<typename Index, typename Scalar, int StorageOrder, bool ConjugateLhs, bool ConjugateRhs> struct general_matrix_vector_product_wrapper
{ staticvoid 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<Scalar,Index,StorageOrder> LhsMapper; typedef internal::const_blas_data_mapper<Scalar,Index,RowMajor> RhsMapper;
const Scalar* a = reinterpret_cast<const Scalar*>(pa); const Scalar* b = reinterpret_cast<const Scalar*>(pb);
Scalar* c = reinterpret_cast<Scalar*>(pc);
Scalar alpha = *reinterpret_cast<const Scalar*>(palpha);
Scalar beta = *reinterpret_cast<const Scalar*>(pbeta);
// check arguments int info = 0; if(OP(*opa)==INVALID) info = 1; elseif(*m<0) info = 2; elseif(*n<0) info = 3; elseif(*lda<std::max(1,*m)) info = 6; elseif(*incb==0) info = 8; elseif(*incc==0) info = 11; if(info) return xerbla_(SCALAR_SUFFIX_UP"GEMV ",&info,6);
const Scalar* a = reinterpret_cast<const Scalar*>(pa);
Scalar* b = reinterpret_cast<Scalar*>(pb);
int info = 0; if(UPLO(*uplo)==INVALID) info = 1; elseif(OP(*opa)==INVALID) info = 2; elseif(DIAG(*diag)==INVALID) info = 3; elseif(*n<0) info = 4; elseif(*lda<std::max(1,*n)) info = 6; elseif(*incb==0) info = 8; if(info) return xerbla_(SCALAR_SUFFIX_UP"TRSV ",&info,6);
const Scalar* a = reinterpret_cast<const Scalar*>(pa);
Scalar* b = reinterpret_cast<Scalar*>(pb);
int info = 0; if(UPLO(*uplo)==INVALID) info = 1; elseif(OP(*opa)==INVALID) info = 2; elseif(DIAG(*diag)==INVALID) info = 3; elseif(*n<0) info = 4; elseif(*lda<std::max(1,*n)) info = 6; elseif(*incb==0) info = 8; if(info) return xerbla_(SCALAR_SUFFIX_UP"TRMV ",&info,6);
/** 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<const Scalar*>(pa); const Scalar* x = reinterpret_cast<const Scalar*>(px);
Scalar* y = reinterpret_cast<Scalar*>(py);
Scalar alpha = *reinterpret_cast<const Scalar*>(palpha);
Scalar beta = *reinterpret_cast<const Scalar*>(pbeta); int coeff_rows = *kl+*ku+1;
int info = 0; if(OP(*trans)==INVALID) info = 1; elseif(*m<0) info = 2; elseif(*n<0) info = 3; elseif(*kl<0) info = 4; elseif(*ku<0) info = 5; elseif(*lda<coeff_rows) info = 8; elseif(*incx==0) info = 10; elseif(*incy==0) info = 13; if(info) return xerbla_(SCALAR_SUFFIX_UP"GBMV ",&info,6);
#if 0 /** TBMV 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 band matrix, with ( k + 1 ) diagonals.
*/ int EIGEN_BLAS_FUNC(tbmv)(char *uplo, char *opa, char *diag, int *n, int *k, RealScalar *pa, int *lda, RealScalar *px, int *incx)
{
Scalar* a = reinterpret_cast<Scalar*>(pa);
Scalar* x = reinterpret_cast<Scalar*>(px); int coeff_rows = *k + 1;
int info = 0; if(UPLO(*uplo)==INVALID) info = 1; elseif(OP(*opa)==INVALID) info = 2; elseif(DIAG(*diag)==INVALID) info = 3; elseif(*n<0) info = 4; elseif(*k<0) info = 5; elseif(*lda<coeff_rows) info = 7; elseif(*incx==0) info = 9; if(info) return xerbla_(SCALAR_SUFFIX_UP"TBMV ",&info,6);
/** DTBSV 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 band matrix, with ( k + 1 ) * diagonals. * * 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(tbsv)(char *uplo, char *op, char *diag, int *n, int *k, RealScalar *pa, int *lda, RealScalar *px, int *incx)
{ typedefvoid (*functype)(int, int, const Scalar *, int, Scalar *); staticconst functype func[16] = { // array index: NOTR | (UP << 2) | (NUNIT << 3)
(internal::band_solve_triangular_selector<int,Upper|0, Scalar,false,Scalar,ColMajor>::run), // array index: TR | (UP << 2) | (NUNIT << 3)
(internal::band_solve_triangular_selector<int,Lower|0, Scalar,false,Scalar,RowMajor>::run), // array index: ADJ | (UP << 2) | (NUNIT << 3)
(internal::band_solve_triangular_selector<int,Lower|0, Scalar,Conj, Scalar,RowMajor>::run),
0, // array index: NOTR | (LO << 2) | (NUNIT << 3)
(internal::band_solve_triangular_selector<int,Lower|0, Scalar,false,Scalar,ColMajor>::run), // array index: TR | (LO << 2) | (NUNIT << 3)
(internal::band_solve_triangular_selector<int,Upper|0, Scalar,false,Scalar,RowMajor>::run), // array index: ADJ | (LO << 2) | (NUNIT << 3)
(internal::band_solve_triangular_selector<int,Upper|0, Scalar,Conj, Scalar,RowMajor>::run),
0, // array index: NOTR | (UP << 2) | (UNIT << 3)
(internal::band_solve_triangular_selector<int,Upper|UnitDiag,Scalar,false,Scalar,ColMajor>::run), // array index: TR | (UP << 2) | (UNIT << 3)
(internal::band_solve_triangular_selector<int,Lower|UnitDiag,Scalar,false,Scalar,RowMajor>::run), // array index: ADJ | (UP << 2) | (UNIT << 3)
(internal::band_solve_triangular_selector<int,Lower|UnitDiag,Scalar,Conj, Scalar,RowMajor>::run),
0, // array index: NOTR | (LO << 2) | (UNIT << 3)
(internal::band_solve_triangular_selector<int,Lower|UnitDiag,Scalar,false,Scalar,ColMajor>::run), // array index: TR | (LO << 2) | (UNIT << 3)
(internal::band_solve_triangular_selector<int,Upper|UnitDiag,Scalar,false,Scalar,RowMajor>::run), // array index: ADJ | (LO << 2) | (UNIT << 3)
(internal::band_solve_triangular_selector<int,Upper|UnitDiag,Scalar,Conj, Scalar,RowMajor>::run),
0,
};
Scalar* a = reinterpret_cast<Scalar*>(pa);
Scalar* x = reinterpret_cast<Scalar*>(px); int coeff_rows = *k+1;
int info = 0; if(UPLO(*uplo)==INVALID) info = 1; elseif(OP(*op)==INVALID) info = 2; elseif(DIAG(*diag)==INVALID) info = 3; elseif(*n<0) info = 4; elseif(*k<0) info = 5; elseif(*lda<coeff_rows) info = 7; elseif(*incx==0) info = 9; if(info) return xerbla_(SCALAR_SUFFIX_UP"TBSV ",&info,6);
/** 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)
{ typedefvoid (*functype)(int, const Scalar*, const Scalar*, Scalar*, Scalar); staticconst functype func[16] = { // array index: NOTR | (UP << 2) | (NUNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Upper|0, Scalar,false,Scalar,false,ColMajor>::run), // array index: TR | (UP << 2) | (NUNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Lower|0, Scalar,false,Scalar,false,RowMajor>::run), // array index: ADJ | (UP << 2) | (NUNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Lower|0, Scalar,Conj, Scalar,false,RowMajor>::run),
0, // array index: NOTR | (LO << 2) | (NUNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Lower|0, Scalar,false,Scalar,false,ColMajor>::run), // array index: TR | (LO << 2) | (NUNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Upper|0, Scalar,false,Scalar,false,RowMajor>::run), // array index: ADJ | (LO << 2) | (NUNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Upper|0, Scalar,Conj, Scalar,false,RowMajor>::run),
0, // array index: NOTR | (UP << 2) | (UNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Upper|UnitDiag,Scalar,false,Scalar,false,ColMajor>::run), // array index: TR | (UP << 2) | (UNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Lower|UnitDiag,Scalar,false,Scalar,false,RowMajor>::run), // array index: ADJ | (UP << 2) | (UNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Lower|UnitDiag,Scalar,Conj, Scalar,false,RowMajor>::run),
0, // array index: NOTR | (LO << 2) | (UNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Lower|UnitDiag,Scalar,false,Scalar,false,ColMajor>::run), // array index: TR | (LO << 2) | (UNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Upper|UnitDiag,Scalar,false,Scalar,false,RowMajor>::run), // array index: ADJ | (LO << 2) | (UNIT << 3)
(internal::packed_triangular_matrix_vector_product<int,Upper|UnitDiag,Scalar,Conj, Scalar,false,RowMajor>::run),
0
};
Scalar* ap = reinterpret_cast<Scalar*>(pap);
Scalar* x = reinterpret_cast<Scalar*>(px);
int info = 0; if(UPLO(*uplo)==INVALID) info = 1; elseif(OP(*opa)==INVALID) info = 2; elseif(DIAG(*diag)==INVALID) info = 3; elseif(*n<0) info = 4; elseif(*incx==0) info = 7; if(info) return xerbla_(SCALAR_SUFFIX_UP"TPMV ",&info,6);
/** 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)
{ typedefvoid (*functype)(int, const Scalar*, Scalar*); staticconst functype func[16] = { // array index: NOTR | (UP << 2) | (NUNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|0, false,ColMajor>::run), // array index: TR | (UP << 2) | (NUNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|0, false,RowMajor>::run), // array index: ADJ | (UP << 2) | (NUNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|0, Conj, RowMajor>::run),
0, // array index: NOTR | (LO << 2) | (NUNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|0, false,ColMajor>::run), // array index: TR | (LO << 2) | (NUNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|0, false,RowMajor>::run), // array index: ADJ | (LO << 2) | (NUNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|0, Conj, RowMajor>::run),
0, // array index: NOTR | (UP << 2) | (UNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|UnitDiag,false,ColMajor>::run), // array index: TR | (UP << 2) | (UNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|UnitDiag,false,RowMajor>::run), // array index: ADJ | (UP << 2) | (UNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|UnitDiag,Conj, RowMajor>::run),
0, // array index: NOTR | (LO << 2) | (UNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|UnitDiag,false,ColMajor>::run), // array index: TR | (LO << 2) | (UNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|UnitDiag,false,RowMajor>::run), // array index: ADJ | (LO << 2) | (UNIT << 3)
(internal::packed_triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|UnitDiag,Conj, RowMajor>::run),
0
};
Scalar* ap = reinterpret_cast<Scalar*>(pap);
Scalar* x = reinterpret_cast<Scalar*>(px);
int info = 0; if(UPLO(*uplo)==INVALID) info = 1; elseif(OP(*opa)==INVALID) info = 2; elseif(DIAG(*diag)==INVALID) info = 3; elseif(*n<0) info = 4; elseif(*incx==0) info = 7; if(info) return xerbla_(SCALAR_SUFFIX_UP"TPSV ",&info,6);
Die Informationen auf dieser Webseite wurden
nach bestem Wissen sorgfältig zusammengestellt. Es wird jedoch weder Vollständigkeit, noch Richtigkeit,
noch Qualität der bereit gestellten Informationen zugesichert.
Bemerkung:
Die farbliche Syntaxdarstellung ist noch experimentell.