?GEMV_

List

sgemv_cblas_sgemv単精度実 数一般行列とベクトルの積
dgemv_cblas_dgemv倍精度実 数一般行列とベクトルの積
cgemv_cblas_cgemv単精度複素数一般行列とベクトルの積
zgemv_cblas_zgemv倍精度複素数一般行列とベクトルの積

概略

一般行列とベクトルの積を計算します。ベクトルが列ベクトルとして解釈される点に注意してください。結果は、渡したベクトルyに格納されます。

計算式

y := alpha * Ax + beta * y

サンプルコード

通常の実数の例
#include <cblas.h>
#include <stdio.h>

int main(void)
{
	int i,j;
	double A[2*3]={1.0,2.0,3.0, 
                       4.0,5.0,6.0};
	double x[3]  ={1.0,2.0,3.0}; //(1 2 3)ベクトル(縦ベクトル)
        double y[2]; //(14 32)となれば良い
        
	cblas_dgemv(CblasRowMajor,//通常のCの行列であればCBlasRowMajorを指定
                    CblasNoTrans, //Aについて転置しない場合 CblasNoTrans 転置する場合  CblasTrans
                    2,            //Aの行数(転置する場合は転置後の)
                    3,            //Aの列数(転置する場合は転置後の)
                    1,            //alphaの値
                    A,            //A
                    3,            //leading dimension(通常はAの列数)
                    x,            //x
                    1,            //incx xのインクリメント幅(通常は1)
                    0,            //betaの値
                    y,            //y
                    1             //incx yのインクリメント幅(通常は1)
                    );
	
	for(j=0;j<2;j++){
		printf("%f ",y[j]);
	}
	printf("\n");
        return 0;
}
Fortranを直呼びする場合
#include <cblas.h>
#include <stdio.h>

int main(void)
{
	int i,j;
	double *A;
	double *x,*y;
	int m=3,n=3;int ONE=1;
	double alpha=1,beta=2;
	
	A=(double*)malloc(3 * 3 * sizeof(double));
	x=(double*)malloc(    3 * sizeof(double));
	y=(double*)malloc(    3 * sizeof(double));
	
	A[0]=1; A[3]=2; A[6]=3;
	A[1]=0; A[4]=2; A[7]=3;
	A[2]=0; A[5]=0; A[8]=3;
	
	x[0]=1;   y[0]=1;
	x[1]=1;   y[1]=1;
	x[2]=1;   y[2]=1;
	
	dgemv_("N",&m,&n,&alpha,A,&m,x,&ONE,&beta,y,&ONE);
	
	for(j=0;j<n;j++){
		printf("%f ",y[j]);
	}
	printf("\n");
        return 0;
}
複素数の例
#include <cblas.h>
#include <stdio.h>

typedef struct
{
  double r;
  double i;
} complex_;

int main(void)
{
  //複素数を構造体無しで扱う場合はaxpyを参照
  complex_ A[2][2]={{{1.0,2.0},{3.0,4.0}},  // 1+2i 3+4i
                    {{5.0,6.0},{7.0,8.0}}}; // 5+6i 7+8i
  complex_ x[2]   = {{3.0,4.0},{6.0,8.0}};  //(3+4i 6+8i)

  complex_ y[2];    //(19+58i -31+142i)となることを確かめる
  complex_ alpha={1,0};
  complex_ beta={0,0};
  int i,j;
        
	cblas_zgemv(CblasRowMajor,//通常のCの行列であればCBlasRowMajorを指定
                    CblasNoTrans, //Aを転置しない場合 CblasNoTrans 転置 CblasTrans 共役転置 CblasConjTrans
                    2,            //Aの行数(転置する場合は転置後の)
                    2,            //Aの列数(転置する場合は転置後の)
                    &alpha,       //alphaの値(構造体へのポインタで指定)
                    A,            //A
                    2,            //leading dimension(通常はAの列数)
                    x,            //x
                    1,            //incx xのインクリメント幅(通常は1)
                    &beta,        //betaの値(構造体へのポインタで指定)
                    y,            //y
                    1             //incx yのインクリメント幅(通常は1)
                    );
	
	for(j=0;j<2;j++){
		printf("%lf %lf i, ",y[j].r,y[j].i);
	}
	printf("\n");
        return 0;
}

引数/戻り値

cblas_dgemv

変数名上書き概要
order enum CBLAS_ORDER C言語では通常はCBlasRowMajor(列方向に並べる)
Fortran式ならCBlasColMajor(行方向に並べる)
trans enum CBLAS_TRANSPOSE 行列の転置を指定
転置しない CblasNoTrans
転置 CblasTrans
共役転置 CblasConjTrans
m int 行列の行数
n int 行列の列数
alpha double スカラーalpha
A double* 行列Aの先頭ポインタ
ldA int Aのleading dimension
CBlasRowMajorを指定したら列数
CBlasColMajor指定なら行数
x double* ベクトルxの先頭ポインタ
incx int Xのインクリメント幅(通常1を指定すれば良い)
beta double スカラーbeta
y double* 上書きベクトルYの先頭ポインタ
incy int* Yのインクリメント幅(通常1を指定すれば良い)
戻り値void

dgemv_

変数名概要
trans char* 行列の転置を指定 ("N"(そのまま),"T"(転置),"C"(共役転置)から選択)
m int* 行列の行数
n int* 行列の列数
alpha double*スカラーalpha
A double*行列Aの先頭ポインタ
ldA int* Aのleading dimension (通常は行数を指定すれば良い)
x double*ベクトルxの先頭ポインタ
incx int* Xのインクリメント幅(通常1を指定すれば良い)
beta double*スカラーbeta
y double*ベクトルYの先頭ポインタ
incy int* Yのインクリメント幅(通常1を指定すれば良い)
戻り値void

プロトタイプ宣言

void sgemv_(char *trans, int *m, int *n,
float *alpha, float *A, int *ldA, float *x, int *incx,
float *beta , float *y, int *incy);

void dgemv_(char *trans, int *m, int *n,
double *alpha, double *A, int *ldA, double *x, int *incx,
double *beta , double *y, int *incy);

void cgemv_(char *trans, int *m, int *n,
complex *alpha, complex *A, int *ldA, complex *x, int *incx,
complex *beta , complex *y, int *incy);

void zgemv_(char *trans, int *m, int *n,
doublecomplex *alpha, doublecomplex *A, int *ldA, doublecomplex *x, int *incx,
doublecomplex *beta , doublecomplex *y, int *incy);

void cblas_sgemv(const enum CBLAS_ORDER order, const enum CBLAS_TRANSPOSE trans, const blasint m, const blasint n,
const float alpha, const float *a, const blasint lda, const float *x, const blasint incx,
const float beta, float *y, const blasint incy);

void cblas_dgemv(const enum CBLAS_ORDER order, const enum CBLAS_TRANSPOSE trans, const blasint m, const blasint n,
const double alpha, const double *a, const blasint lda, const double *x, const blasint incx,
const double beta, double *y, const blasint incy);

void cblas_cgemv(const enum CBLAS_ORDER order, const enum CBLAS_TRANSPOSE trans, const blasint m, const blasint n,
const void *alpha, const void *a, const blasint lda, const void *x, const blasint incx,
const void *beta, void *y, const blasint incy);

void cblas_zgemv(const enum CBLAS_ORDER order, const enum CBLAS_TRANSPOSE trans, const blasint m, const blasint n,
const void *alpha, const void *a, const blasint lda, const void *x, const blasint incx,
const void *beta, void *y, const blasint incy);