[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

13 Suporte a BLAS

O Basic Álgebra Linear Subprograms (BLAS) define um conjunto de operações fundamentais sobre vetores e matrizes que podem ser usados para criar funcionalidades otimizadas de alto nível para serem usadas em álgebra linear.

A biblioteca fornece uma camada de baixo nível que corresponde diretamente à linguagem C padrão usada pela BLAS padrão, referida aqui como “CBLAS”, e uma interface de alto nível para operações sobre vetores e matrizes da GSL. Usuários que estão interessados em operações simples sobre vetores e matrizes em quanto objetos GSL devem usar a camada de alto nível descrita nesse capítulo. As funções são declaradas no arquivo ‘gsl_blas.h’ e deve satisfazer as necessidades da maioria dos usuários.

Note que matrizes GSL estão implementadas usando armazenagem densa de forma que a interface somente inclui as funções correspondentes de armazenagem densa BLAS. A funcionalidade plena da BLAS para armazenagem de matrizes em banda e armazenagem de matrizes em pacote está disponível através da interface CBLAS de baixo nível. Similarmente, vetores da GSL estão restritos a saltos positivos, ao passo que a interface CBLAS de baixo nível suporta saltos negativos como especificado no padrão BLAS.(34)

A interface para a camada gsl_cblas é especificado no arquivo ‘gsl_cblas.h’. Essa interface corresponde ao padrão de fóruns técnicos da BLAS para a interface C para implementações de BLAS herdadas. Usuários que possuem acesso a outras implementações CBLAS que seguem as regras padronizadas podem usar ‘gsl_cblas.h’ em lugar de versões fornecidas pela biblioteca. Note que usuários que possuem somente uma biblioteca BLAS de Fortran podem usar um pacote CBLAS que obedece às especificações para converter seus arquivos em Fortran para uma biblioteca CBLAS. Uma referência a pacote CBLAS para implemetações Fortran herdadas existe como parte da CBLAS padrão e pode ser obtida na Netlib (35). O conjunto completo de funções da CBLAS é listada em um apêndice (veja seção Biblioteca CBLAS da GSL).

Existem três níveis de operações BLAS,

Nível 1

Operações com vetores, e.g. y = \alpha x + y

Nível 2

Operações matriz-vetor, e.g. y = \alpha A x + \beta y

Nível 3

Operações matriz-matriz, e.g. C = \alpha A B + C

Cada rotina tem um nome que especifica a operação, o tipo de matrizes involvidas e sua precisão. Algumas das operações mais comuns e seus nomes são fornecidas abaixo,

DOT

produto escalar, x^T y

AXPY

soma de vetores, \alpha x + y

MV

produto matriz-vetor, A x

SV

resolve matriz-vetor, inv(A) x

MM

produto matriz-matriz, A B

SM

resolve matriz-matriz, inv(A) B

Os tipos de matrizes são,

GE

geral

GB

banda geral

SY

simétrica

SB

banda simétrica

SP

pacote simétrica

HE

de Hermite

HB

banda de Hermite

HP

pacote de Hermite

TR

triangular

TB

banda triangular

TP

pacote triangular

Cada operação está definida para quatro precisões,

S

real de precisão simples

D

real de precisão dupla

C

complexo de precisão simples

Z

complexo de precisão dupla

Dessa forma, por exemplo, o nome SGEMM significa “single-precision general matrix-matrix multiply” e ZGEMM significa “double-precision complex matrix-matrix multiply”.

Note que argumentos vetores e matrizes para funções BLAS não devem ter nomes alternativos (alias), já que os resultados são indefinidos quando os respectivos vetores estáticos sobrescrevem (veja seção Nomes alternativos para vetores estáticos).


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

13.1 Interface GSL BLAS

GSL fornece objetos densos em forma de vetores e em forma de matrizes, baseado sobre tipos internos relevantes. A biblioteca fornece uma interface para operações na BLAS que aplicam esses objetos. A interface para essa funcionalidade é fornecida no arquivo ‘gsl_blas.h’.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

13.1.1 Nível 1

Function: int gsl_blas_sdsdot (float alpha, const gsl_vector_float * x, const gsl_vector_float * y, float * result)

Essa função calcula a soma \alpha + x^T y para os vetores x e y, retornando o resultado em result.

Function: int gsl_blas_sdot (const gsl_vector_float * x, const gsl_vector_float * y, float * result)
Function: int gsl_blas_dsdot (const gsl_vector_float * x, const gsl_vector_float * y, double * result)
Function: int gsl_blas_ddot (const gsl_vector * x, const gsl_vector * y, double * result)

Essas funções calculam o produto escalar x^T y para os vetores x e y, retornando o resultado em result.

Function: int gsl_blas_cdotu (const gsl_vector_complex_float * x, const gsl_vector_complex_float * y, gsl_complex_float * dotu)
Function: int gsl_blas_zdotu (const gsl_vector_complex * x, const gsl_vector_complex * y, gsl_complex * dotu)

Essas funções calculam o produto escalar complexo x^T y para os vetores x e y, retornando o resultado em dotu

Function: int gsl_blas_cdotc (const gsl_vector_complex_float * x, const gsl_vector_complex_float * y, gsl_complex_float * dotc)
Function: int gsl_blas_zdotc (const gsl_vector_complex * x, const gsl_vector_complex * y, gsl_complex * dotc)

Essas funções calculam o produto escalar conjugado complexo x^H y para os vetores x e y, retornando o resutado em dotc

Function: float gsl_blas_snrm2 (const gsl_vector_float * x)
Function: double gsl_blas_dnrm2 (const gsl_vector * x)

Essas funções calculam o módulo ou norma Euclideana ||x||_2 = \sqrt {\sum x_i^2} do vetor x.

Function: float gsl_blas_scnrm2 (const gsl_vector_complex_float * x)
Function: double gsl_blas_dznrm2 (const gsl_vector_complex * x)

Essas funções calculam o módulo ou norma Euclideana do vetor complexo x,

||x||2 =   ⎛



(Re(xi)2 + Im(xi)2)
 
.

Function: float gsl_blas_sasum (const gsl_vector_float * x)
Function: double gsl_blas_dasum (const gsl_vector * x)

Essas funções calculam a soma absoluta \sum |x_i| dos elementos do vetor x.

Function: float gsl_blas_scasum (const gsl_vector_complex_float * x)
Function: double gsl_blas_dzasum (const gsl_vector_complex * x)

Essas funções calculam a soma dos valores absolutos das partes real e imaginárias do vetor complexo x, \sum |\Re(x_i)| + |\Im(x_i)|.

Function: CBLAS_INDEX_t gsl_blas_isamax (const gsl_vector_float * x)
Function: CBLAS_INDEX_t gsl_blas_idamax (const gsl_vector * x)
Function: CBLAS_INDEX_t gsl_blas_icamax (const gsl_vector_complex_float * x)
Function: CBLAS_INDEX_t gsl_blas_izamax (const gsl_vector_complex * x)

Essas funções retornam o índice do maior elemento do vetor x. O maior elemento é determinado por seu valor absoluto para vetores reais e pela soma dos valores absolutos das partes reais e imaginárias |\Re(x_i)| + |\Im(x_i)| para vetores complexos. Se o maior valor ocorrer repetidas vezes então o índice da primeira ocorrência é retornado.

Function: int gsl_blas_sswap (gsl_vector_float * x, gsl_vector_float * y)
Function: int gsl_blas_dswap (gsl_vector * x, gsl_vector * y)
Function: int gsl_blas_cswap (gsl_vector_complex_float * x, gsl_vector_complex_float * y)
Function: int gsl_blas_zswap (gsl_vector_complex * x, gsl_vector_complex * y)

Essas funções trocam os elementos dos vetores x e y.

Function: int gsl_blas_scopy (const gsl_vector_float * x, gsl_vector_float * y)
Function: int gsl_blas_dcopy (const gsl_vector * x, gsl_vector * y)
Function: int gsl_blas_ccopy (const gsl_vector_complex_float * x, gsl_vector_complex_float * y)
Function: int gsl_blas_zcopy (const gsl_vector_complex * x, gsl_vector_complex * y)

Essas funções copiam os elementos do vetor x para o vetor y.

Function: int gsl_blas_saxpy (float alpha, const gsl_vector_float * x, gsl_vector_float * y)
Function: int gsl_blas_daxpy (double alpha, const gsl_vector * x, gsl_vector * y)
Function: int gsl_blas_caxpy (const gsl_complex_float alpha, const gsl_vector_complex_float * x, gsl_vector_complex_float * y)
Function: int gsl_blas_zaxpy (const gsl_complex alpha, const gsl_vector_complex * x, gsl_vector_complex * y)

Essas funções calculam a soma y = \alpha x + y para os vetores x e y.

Function: void gsl_blas_sscal (float alpha, gsl_vector_float * x)
Function: void gsl_blas_dscal (double alpha, gsl_vector * x)
Function: void gsl_blas_cscal (const gsl_complex_float alpha, gsl_vector_complex_float * x)
Function: void gsl_blas_zscal (const gsl_complex alpha, gsl_vector_complex * x)
Function: void gsl_blas_csscal (float alpha, gsl_vector_complex_float * x)
Function: void gsl_blas_zdscal (double alpha, gsl_vector_complex * x)

Essas funções alteram o vetor x através de um fator multiplicativo alpha.

Function: int gsl_blas_srotg (float a[], float b[], float c[], float s[])
Function: int gsl_blas_drotg (double a[], double b[], double c[], double s[])

Essas funções calculam uma rotação de Givens (c,s) que zera o vetor (a,b),




c
s
−s
c






a
b



=


r′
0



As variáveis a e b são sobrescritas pela rotina.

Function: int gsl_blas_srot (gsl_vector_float * x, gsl_vector_float * y, float c, float s)
Function: int gsl_blas_drot (gsl_vector * x, gsl_vector * y, const double c, const double s)

Essas funções aplicam uma rotação de Givens (x', y') = (c x + s y, -s x + c y) aos vetores x, y.

Function: int gsl_blas_srotmg (float d1[], float d2[], float b1[], float b2, float P[])
Function: int gsl_blas_drotmg (double d1[], double d2[], double b1[], double b2, double P[])

Essas funções calculam uma transformação de Givens modificada. A transformação de Givens modificada é definida na especificação da BLAS Nível 1 original, fornecida nas referências.

Function: int gsl_blas_srotm (gsl_vector_float * x, gsl_vector_float * y, const float P[])
Function: int gsl_blas_drotm (gsl_vector * x, gsl_vector * y, const double P[])

Essas funções aplicam uma transformação de Givens.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

13.1.2 Nível 2

Function: int gsl_blas_sgemv (CBLAS_TRANSPOSE_t TransA, float alpha, const gsl_matrix_float * A, const gsl_vector_float * x, float beta, gsl_vector_float * y)
Function: int gsl_blas_dgemv (CBLAS_TRANSPOSE_t TransA, double alpha, const gsl_matrix * A, const gsl_vector * x, double beta, gsl_vector * y)
Function: int gsl_blas_cgemv (CBLAS_TRANSPOSE_t TransA, const gsl_complex_float alpha, const gsl_matrix_complex_float * A, const gsl_vector_complex_float * x, const gsl_complex_float beta, gsl_vector_complex_float * y)
Function: int gsl_blas_zgemv (CBLAS_TRANSPOSE_t TransA, const gsl_complex alpha, const gsl_matrix_complex * A, const gsl_vector_complex * x, const gsl_complex beta, gsl_vector_complex * y)

Essas funções calculam a soma e o produto matriz-vetor y = \alpha op(A) x + \beta y, onde op(A) = A, A^T, A^H para TransA = CblasNoTrans, CblasTrans, CblasConjTrans.

Function: int gsl_blas_strmv (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_matrix_float * A, gsl_vector_float * x)
Function: int gsl_blas_dtrmv (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_matrix * A, gsl_vector * x)
Function: int gsl_blas_ctrmv (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_matrix_complex_float * A, gsl_vector_complex_float * x)
Function: int gsl_blas_ztrmv (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_matrix_complex * A, gsl_vector_complex * x)

Essas funções calculam o produto matriz-vetor x = op(A) x para a matriz triangular A, onde op(A) = A, A^T, A^H para TransA = CblasNoTrans, CblasTrans, CblasConjTrans. Quando Uplo for CblasUpper então o triângulo alto de A é usado, e quando Uplo for CblasLower então o triângulo baixo de A é usado. Se Diag for CblasNonUnit então a diagonal da matriz é usada, mas se Diag for CblasUnit então os elementos da diagonal da matriz A são tomados como unitários e não são referenciados.

Function: int gsl_blas_strsv (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_matrix_float * A, gsl_vector_float * x)
Function: int gsl_blas_dtrsv (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_matrix * A, gsl_vector * x)
Function: int gsl_blas_ctrsv (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_matrix_complex_float * A, gsl_vector_complex_float * x)
Function: int gsl_blas_ztrsv (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_matrix_complex * A, gsl_vector_complex * x)

Essas funções calculam inv(op(A)) x para x, onde op(A) = A, A^T, A^H para TransA = CblasNoTrans, CblasTrans, CblasConjTrans. Quando Uplo for CblasUpper então o triângulo alto de A é usado, e quando Uplo for CblasLower então o triângulo baixo de A é usado. Se Diag for CblasNonUnit então a diagonal da matriz é usada, mas se Diag for CblasUnit então os elementos da diagonal da matriz A são tomados como sendo unitários e não são referenciados.

Function: int gsl_blas_ssymv (CBLAS_UPLO_t Uplo, float alpha, const gsl_matrix_float * A, const gsl_vector_float * x, float beta, gsl_vector_float * y)
Function: int gsl_blas_dsymv (CBLAS_UPLO_t Uplo, double alpha, const gsl_matrix * A, const gsl_vector * x, double beta, gsl_vector * y)

Essas funções calculam a soma e o produto matriz-vetor y = \alpha A x + \beta y para a matriz simétrica A. Uma vez que a matriz A é simétrica somente sua meia parte alta ou baixa precisa ser armazenada. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de A são usados, e quando Uplo for CblasLower então o triângulo baixo e a diagonal de A são usados.

Function: int gsl_blas_chemv (CBLAS_UPLO_t Uplo, const gsl_complex_float alpha, const gsl_matrix_complex_float * A, const gsl_vector_complex_float * x, const gsl_complex_float beta, gsl_vector_complex_float * y)
Function: int gsl_blas_zhemv (CBLAS_UPLO_t Uplo, const gsl_complex alpha, const gsl_matrix_complex * A, const gsl_vector_complex * x, const gsl_complex beta, gsl_vector_complex * y)

Essas funções calculam o produto e a soma matriz-vetor y = \alpha A x + \beta y para a matriz de Hermite A. Uma vez que a matriz A é de Hermite somente sua metade alta ou baixa precisa ser armazenada. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de A são usados, e quando Uplo for CblasLower então o triângulo baixo e a diagonal de A são usados. Os elementos imaginários da diagonal são automaticamente assumidos serem zero e não são referenciados.

Function: int gsl_blas_sger (float alpha, const gsl_vector_float * x, const gsl_vector_float * y, gsl_matrix_float * A)
Function: int gsl_blas_dger (double alpha, const gsl_vector * x, const gsl_vector * y, gsl_matrix * A)
Function: int gsl_blas_cgeru (const gsl_complex_float alpha, const gsl_vector_complex_float * x, const gsl_vector_complex_float * y, gsl_matrix_complex_float * A)
Function: int gsl_blas_zgeru (const gsl_complex alpha, const gsl_vector_complex * x, const gsl_vector_complex * y, gsl_matrix_complex * A)

Essas funções calculam a atualização rank-1 A = \alpha x y^T + A da matriz A.

Function: int gsl_blas_cgerc (const gsl_complex_float alpha, const gsl_vector_complex_float * x, const gsl_vector_complex_float * y, gsl_matrix_complex_float * A)
Function: int gsl_blas_zgerc (const gsl_complex alpha, const gsl_vector_complex * x, const gsl_vector_complex * y, gsl_matrix_complex * A)

Essas funções calculam o conjugado da atualização rank-1 A = \alpha x y^H + A da matriz A.

Function: int gsl_blas_ssyr (CBLAS_UPLO_t Uplo, float alpha, const gsl_vector_float * x, gsl_matrix_float * A)
Function: int gsl_blas_dsyr (CBLAS_UPLO_t Uplo, double alpha, const gsl_vector * x, gsl_matrix * A)

Essas funções calculam a atualização rank-1 simétrica A = \alpha x x^T + A da matriz simétrica A. Uma vez que a matriz A é simétrica somente sua metade alta ou baixa precisa ser armazenada. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de A são usados, e quando Uplo for CblasLower então o triângulo baixo e a diagonal de A são usados.

Function: int gsl_blas_cher (CBLAS_UPLO_t Uplo, float alpha, const gsl_vector_complex_float * x, gsl_matrix_complex_float * A)
Function: int gsl_blas_zher (CBLAS_UPLO_t Uplo, double alpha, const gsl_vector_complex * x, gsl_matrix_complex * A)

Essas funções calculam a atualização rank-1 de Hermite A = \alpha x x^H + A da matriz de Hermite A. Uma vez que a matriz A é de Hermite somente sua metade alta ou baixa precisa ser armazenada. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de A são usados, e quando Uplo for CblasLower então o triângulo baixo e a diagonal de A são usados. Os elementos imaginários da diagonal são automaticamente ajustados para zero.

Function: int gsl_blas_ssyr2 (CBLAS_UPLO_t Uplo, float alpha, const gsl_vector_float * x, const gsl_vector_float * y, gsl_matrix_float * A)
Function: int gsl_blas_dsyr2 (CBLAS_UPLO_t Uplo, double alpha, const gsl_vector * x, const gsl_vector * y, gsl_matrix * A)

Essas funções calculam a atualização rank-2 simétrica A = \alpha x y^T + \alpha y x^T + A da matriz simétrica A. Uma vez que a matriz A é simétrica somente sua metade alta ou baixa precisa ser armazenada. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de A são usados, e quando Uplo for CblasLower enttão o triângulo baixo e a diagonal de A são usados.

Function: int gsl_blas_cher2 (CBLAS_UPLO_t Uplo, const gsl_complex_float alpha, const gsl_vector_complex_float * x, const gsl_vector_complex_float * y, gsl_matrix_complex_float * A)
Function: int gsl_blas_zher2 (CBLAS_UPLO_t Uplo, const gsl_complex alpha, const gsl_vector_complex * x, const gsl_vector_complex * y, gsl_matrix_complex * A)

Essas funções calculam a atualização rank-2 de Hermite A = \alpha x y^H + \alpha^* y x^H + A da matriz de Hermite A. Uma vez que a matriz A é de Hermite somente sua metade alta ou baixa precisa ser armazenada. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de A são usados, e quando Uplo for CblasLower então o triângulo baixo e a diagonal de A são usados. Os elementos imaginários da diagona são automaticamente ajustados para zero.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

13.1.3 Nível 3

Function: int gsl_blas_sgemm (CBLAS_TRANSPOSE_t TransA, CBLAS_TRANSPOSE_t TransB, float alpha, const gsl_matrix_float * A, const gsl_matrix_float * B, float beta, gsl_matrix_float * C)
Function: int gsl_blas_dgemm (CBLAS_TRANSPOSE_t TransA, CBLAS_TRANSPOSE_t TransB, double alpha, const gsl_matrix * A, const gsl_matrix * B, double beta, gsl_matrix * C)
Function: int gsl_blas_cgemm (CBLAS_TRANSPOSE_t TransA, CBLAS_TRANSPOSE_t TransB, const gsl_complex_float alpha, const gsl_matrix_complex_float * A, const gsl_matrix_complex_float * B, const gsl_complex_float beta, gsl_matrix_complex_float * C)
Function: int gsl_blas_zgemm (CBLAS_TRANSPOSE_t TransA, CBLAS_TRANSPOSE_t TransB, const gsl_complex alpha, const gsl_matrix_complex * A, const gsl_matrix_complex * B, const gsl_complex beta, gsl_matrix_complex * C)

Essas funções calculam o produto matriz-matriz e a soma C = \alpha op(A) op(B) + \beta C onde op(A) = A, A^T, A^H para TransA = CblasNoTrans, CblasTrans, CblasConjTrans e similarmente para o parâmetro TransB.

Function: int gsl_blas_ssymm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, float alpha, const gsl_matrix_float * A, const gsl_matrix_float * B, float beta, gsl_matrix_float * C)
Function: int gsl_blas_dsymm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, double alpha, const gsl_matrix * A, const gsl_matrix * B, double beta, gsl_matrix * C)
Function: int gsl_blas_csymm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, const gsl_complex_float alpha, const gsl_matrix_complex_float * A, const gsl_matrix_complex_float * B, const gsl_complex_float beta, gsl_matrix_complex_float * C)
Function: int gsl_blas_zsymm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, const gsl_complex alpha, const gsl_matrix_complex * A, const gsl_matrix_complex * B, const gsl_complex beta, gsl_matrix_complex * C)

Essas funções calculam o produto matriz-matriz e a soma C = \alpha A B + \beta C para Side sendo CblasLeft e C = \alpha B A + \beta C para Side sendo CblasRight, onde a matriz A é simétrica. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de A são usados, e quando Uplo for CblasLower então o triângulo baixo e a diagonal de A são usados.

Function: int gsl_blas_chemm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, const gsl_complex_float alpha, const gsl_matrix_complex_float * A, const gsl_matrix_complex_float * B, const gsl_complex_float beta, gsl_matrix_complex_float * C)
Function: int gsl_blas_zhemm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, const gsl_complex alpha, const gsl_matrix_complex * A, const gsl_matrix_complex * B, const gsl_complex beta, gsl_matrix_complex * C)

Essas funções calculam o produto matriz-matriz e a soma C = \alpha A B + \beta C para Side sendo CblasLeft e C = \alpha B A + \beta C para Side sendo CblasRight, onde a matriz A é de Hermite. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de A são usados, e quando Uplo for CblasLower então o triângulo baixo e a diagonal de A são usados. Os elementos imaginários da diagonal são automaticamente ajustados para zero.

Function: int gsl_blas_strmm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, float alpha, const gsl_matrix_float * A, gsl_matrix_float * B)
Function: int gsl_blas_dtrmm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, double alpha, const gsl_matrix * A, gsl_matrix * B)
Function: int gsl_blas_ctrmm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_complex_float alpha, const gsl_matrix_complex_float * A, gsl_matrix_complex_float * B)
Function: int gsl_blas_ztrmm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_complex alpha, const gsl_matrix_complex * A, gsl_matrix_complex * B)

Essas funções calculam o produto matriz-matriz B = \alpha op(A) B para Side sendo CblasLeft e B = \alpha B op(A) para Side sendo CblasRight. A matriz A é triangular e op(A) = A, A^T, A^H para TransA = CblasNoTrans, CblasTrans, CblasConjTrans. Quando Uplo for CblasUpper então o triângulo alto de A é usado, e quando Uplo for CblasLower então o triângulo baixo de A é usado. Se Diag for CblasNonUnit então a diagonal de A é usada, mas se Diag for CblasUnit então os elementos da diagonal da matriz A são tomados como unitários e não são referenciados.

Function: int gsl_blas_strsm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, float alpha, const gsl_matrix_float * A, gsl_matrix_float * B)
Function: int gsl_blas_dtrsm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, double alpha, const gsl_matrix * A, gsl_matrix * B)
Function: int gsl_blas_ctrsm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_complex_float alpha, const gsl_matrix_complex_float * A, gsl_matrix_complex_float * B)
Function: int gsl_blas_ztrsm (CBLAS_SIDE_t Side, CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t TransA, CBLAS_DIAG_t Diag, const gsl_complex alpha, const gsl_matrix_complex * A, gsl_matrix_complex * B)

Essas funções calculam produto matriz inversa matriz B = \alpha op(inv(A))B para Side sendo CblasLeft e B = \alpha B op(inv(A)) para Side sendo CblasRight. A matriz A é triangular e op(A) = A, A^T, A^H para TransA = CblasNoTrans, CblasTrans, CblasConjTrans. Quando Uplo for CblasUpper então o triângulo alto de A é usado, e quando Uplo for CblasLower então o triângulo baixo de A é usado. Se Diag for CblasNonUnit então a diagonal de A é usada, mas se Diag for CblasUnit então os elementos da diagonal da matriz A são tomados como unitários e não são referenciados.

Function: int gsl_blas_ssyrk (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, float alpha, const gsl_matrix_float * A, float beta, gsl_matrix_float * C)
Function: int gsl_blas_dsyrk (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, double alpha, const gsl_matrix * A, double beta, gsl_matrix * C)
Function: int gsl_blas_csyrk (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, const gsl_complex_float alpha, const gsl_matrix_complex_float * A, const gsl_complex_float beta, gsl_matrix_complex_float * C)
Function: int gsl_blas_zsyrk (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, const gsl_complex alpha, const gsl_matrix_complex * A, const gsl_complex beta, gsl_matrix_complex * C)

Essas funções calcula a atualização rank-k da matriz simétrica C, C = \alpha A A^T + \beta C quando Trans for CblasNoTrans e C = \alpha A^T A + \beta C quando Trans for CblasTrans. Uma vez que a matriz C é simétrica somente sua metade alta ou baixa precisa ser armazenada. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de C são usados, e quando Uplo for CblasLower então o triângulo baixo e a diagonal de C são usados.

Function: int gsl_blas_cherk (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, float alpha, const gsl_matrix_complex_float * A, float beta, gsl_matrix_complex_float * C)
Function: int gsl_blas_zherk (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, double alpha, const gsl_matrix_complex * A, double beta, gsl_matrix_complex * C)

Essas funções calculam uma atualização rank-k da matriz de Hermite C, C = \alpha A A^H + \beta C quando Trans for CblasNoTrans e C = \alpha A^H A + \beta C quando Trans for CblasConjTrans. Uma vez que a matriz C é de Hermite somente sua metade alta ou baixa precisa ser armazenada. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de C são usados, e quando Uplo for CblasLower então o triângulo baixo e a diagonal de C são usados. Os elementos imaginários da diagonal são automaticamente ajustados para zero.

Function: int gsl_blas_ssyr2k (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, float alpha, const gsl_matrix_float * A, const gsl_matrix_float * B, float beta, gsl_matrix_float * C)
Function: int gsl_blas_dsyr2k (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, double alpha, const gsl_matrix * A, const gsl_matrix * B, double beta, gsl_matrix * C)
Function: int gsl_blas_csyr2k (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, const gsl_complex_float alpha, const gsl_matrix_complex_float * A, const gsl_matrix_complex_float * B, const gsl_complex_float beta, gsl_matrix_complex_float * C)
Function: int gsl_blas_zsyr2k (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, const gsl_complex alpha, const gsl_matrix_complex * A, const gsl_matrix_complex * B, const gsl_complex beta, gsl_matrix_complex * C)

Essas funções calculam uma atualização rank-2k da matriz simétrica C, C = \alpha A B^T + \alpha B A^T + \beta C quando Trans for CblasNoTrans e C = \alpha A^T B + \alpha B^T A + \beta C quando Trans for CblasTrans. Uma vez que a matriz C é simétrica somente sua metade alta ou baixa precisa ser armazenada. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de C são usados, e quando Uplo for CblasLower então o trângulo baixo e a diagonal de C são usados.

Function: int gsl_blas_cher2k (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, const gsl_complex_float alpha, const gsl_matrix_complex_float * A, const gsl_matrix_complex_float * B, float beta, gsl_matrix_complex_float * C)
Function: int gsl_blas_zher2k (CBLAS_UPLO_t Uplo, CBLAS_TRANSPOSE_t Trans, const gsl_complex alpha, const gsl_matrix_complex * A, const gsl_matrix_complex * B, double beta, gsl_matrix_complex * C)

Essas funções calculam uma atualização rank-2k da matriz de Hermite C, C = \alpha A B^H + \alpha^* B A^H + \beta C quando Trans for CblasNoTrans e C = \alpha A^H B + \alpha^* B^H A + \beta C quando Trans for CblasConjTrans. Uma vez que a matriz C é de Hermite somente sua metade alta ou baixa precisa ser armazenada. Quando Uplo for CblasUpper então o triângulo alto e a diagonal de C são usados, e quando Uplo for CblasLower então o triângulo baixo e a diagonal de C são usados. Os elementos imaginários da diagonal são automaticamente ajustados para zero.


[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

13.2 Exemplos

O seguinte programa calcula o produto de duas matrizes usando a função DGEMM da BLAS Nível-3,




0.11
0.12
0.13
0.21
0.22
0.23







1011
1012
1021
1022
1031
1031




=


367.76
368.12
674.06
674.72



As matrizes são armazenadas na ordem que privilegia a linha (row major order) , conforme a convenção na linguagem C para vetores estáticos.

#include <stdio.h>
#include <gsl/gsl_blas.h>

int
main (void)
{
  double a[] = { 0.11, 0.12, 0.13,
                 0.21, 0.22, 0.23 };

  double b[] = { 1011, 1012,
                 1021, 1022,
                 1031, 1032 };

  double c[] = { 0.00, 0.00,
                 0.00, 0.00 };

  gsl_matrix_view A = gsl_matrix_view_array(a, 2, 3);
  gsl_matrix_view B = gsl_matrix_view_array(b, 3, 2);
  gsl_matrix_view C = gsl_matrix_view_array(c, 2, 2);

  /* Compute C = A B */

  gsl_blas_dgemm (CblasNoTrans, CblasNoTrans,
                  1.0, &A.matrix, &B.matrix,
                  0.0, &C.matrix);

  printf ("[ %g, %g\n", c[0], c[1]);
  printf ("  %g, %g ]\n", c[2], c[3]);

  return 0;  
}

Aqui está a saída do programa,

$ ./a.out
[ 367.76, 368.12
  674.06, 674.72 ]

[ << ] [ < ] [ Acima ] [ > ] [ >> ]         [Topo] [Conteúdo] [Índice] [ ? ]

13.3 Referências e Leituras Adicionais

Informação sobre as especificações da BLAS, incluindo ambas as especificações da interface herdada e as especificações da interface atualizada, está disponível online no sítio da BLAS e no sítio do BLAS Technical Forum.

Os seguintes artigos possuem especificações para os níveis da BLAS Nível 1, Nível 2 e Nível 3.

Versões no formato postscript dos últimos dois artigos estão disponíveis em http://www.netlib.org/blas/. Um pacote CBLAS para bibliotecas Fortran BLAS está disponível na mesma localização.


[ << ] [ >> ]           [Topo] [Conteúdo] [Índice] [ ? ]

Esse documento foi gerado em 23 de Julho de 2013 usando texi2html 5.0.