|
|
| |
Math::GSL::BLAS(3) |
User Contributed Perl Documentation |
Math::GSL::BLAS(3) |
Math::GSL::BLAS - Basic Linear Algebra Subprograms
use Math::GSL::BLAS qw/:all/;
use Math::GSL::Matrix qw/:all/;
# matrix-matrix product of double numbers
my $A = Math::GSL::Matrix->new(2,2);
$A->set_row(0, [1, 4]);
->set_row(1, [3, 2]);
my $B = Math::GSL::Matrix->new(2,2);
$B->set_row(0, [2, 1]);
->set_row(1, [5,3]);
my $C = Math::GSL::Matrix->new(2,2);
gsl_matrix_set_zero($C->raw);
gsl_blas_dgemm($CblasNoTrans, $CblasNoTrans, 1, $A->raw, $B->raw, 1, $C->raw);
my @got = $C->row(0)->as_list;
print "The resulting matrix is: \n[";
print "$got[0] $got[1]\n";
@got = $C->row(1)->as_list;
print "$got[0] $got[1] ]\n";
# compute the scalar product of two vectors :
use Math::GSL::Vector qw/:all/;
use Math::GSL::CBLAS qw/:all/;
my $vec1 = Math::GSL::Vector->new([1,2,3,4,5]);
my $vec2 = Math::GSL::Vector->new([5,4,3,2,1]);
my ($status, $result) = gsl_blas_ddot($vec1->raw, $vec2->raw);
if($status == 0) {
print "The function has succeeded. \n";
}
print "The result of the vector multiplication is $result.\n";
The functions of this module are divised into 3 levels:
- "gsl_blas_sdsdot"
- "gsl_blas_dsdot"
- "gsl_blas_sdot"
- "gsl_blas_ddot($x, $y)"
- This function computes the scalar product x^T y for the vectors
$x and $y. The function
returns two values, the first is 0 if the operation suceeded, 1 otherwise
and the second value is the result of the computation.
- "gsl_blas_cdotu"
- "gsl_blas_cdotc"
- "gsl_blas_zdotu($x, $y, $dotu)"
- This function computes the complex scalar product x^T y for the complex
vectors $x and $y,
returning the result in the complex number $dotu.
The function returns 0 if the operation suceeded, 1 otherwise.
- "gsl_blas_zdotc($x, $y, $dotc)"
- This function computes the complex conjugate scalar product x^H y for the
complex vectors $x and $y,
returning the result in the complex number $dotc.
The function returns 0 if the operation suceeded, 1 otherwise.
- "gsl_blas_snrm2" =item "gsl_blas_sasum"
- "gsl_blas_dnrm2($x)"
- This function computes the Euclidean norm
||x||_2 = \sqrt {\sum x_i^2}
of the vector $x.
- "gsl_blas_dasum($x)"
- This function computes the absolute sum \sum |x_i| of the elements of the
vector $x.
- "gsl_blas_scnrm2"
- "gsl_blas_scasum"
- "gsl_blas_dznrm2($x)"
- This function computes the Euclidean norm of the complex vector
$x, ||x||_2 = \sqrt {\sum (\Re(x_i)^2 +
\Im(x_i)^2)}.
- "gsl_blas_dzasum($x)"
- This function computes the sum of the magnitudes of the real and imaginary
parts of the complex vector $x, \sum |\Re(x_i)| +
|\Im(x_i)|.
- "gsl_blas_isamax"
- "gsl_blas_idamax"
- "gsl_blas_icamax"
- "gsl_blas_izamax "
- "gsl_blas_sswap"
- "gsl_blas_scopy"
- "gsl_blas_saxpy"
- "gsl_blas_dswap($x, $y)"
- This function exchanges the elements of the vectors
$x and $y. The function
returns 0 if the operation suceeded, 1 otherwise.
- "gsl_blas_dcopy($x, $y)"
- This function copies the elements of the vector $x
into the vector $y. The function returns 0 if the
operation suceeded, 1 otherwise.
- "gsl_blas_daxpy($alpha, $x, $y)"
- These functions compute the sum $y =
$alpha * $x +
$y for the vectors $x and
$y.
- "gsl_blas_cswap"
- "gsl_blas_ccopy "
- "gsl_blas_caxpy"
- "gsl_blas_zswap"
- "gsl_blas_zcopy"
- "gsl_blas_zaxpy "
- "gsl_blas_srotg"
- "gsl_blas_srotmg"
- "gsl_blas_srot"
- "gsl_blas_srotm "
- "gsl_blas_drotg"
- "gsl_blas_drotmg"
- "gsl_blas_drot($x, $y, $c, $s)"
- This function applies a Givens rotation (x', y') = (c x + s y, -s x + c y)
to the vectors $x,
$y.
- "gsl_blas_drotm "
- "gsl_blas_sscal"
- "gsl_blas_dscal($alpha, $x)"
- This function rescales the vector $x by the
multiplicative factor $alpha.
- "gsl_blas_cscal"
- "gsl_blas_zscal "
- "gsl_blas_csscal"
- "gsl_blas_zdscal"
- "gsl_blas_sgemv"
- "gsl_blas_strmv "
- "gsl_blas_strsv"
- "gsl_blas_dgemv($TransA, $alpha, $A, $x, $beta, $y)" - This
function computes the matrix-vector product and sum y = \alpha op(A) x +
\beta y, where op(A) = A, A^T, A^H for $TransA = $CblasNoTrans, $CblasTrans,
$CblasConjTrans (constant values coming from the CBLAS module). $A is a
matrix and $x and $y are vectors. The function returns 0 if the operation
suceeded, 1 otherwise.
- "gsl_blas_dtrmv($Uplo, $TransA, $Diag, $A, $x)" - This function
computes the matrix-vector product x = op(A) x for the triangular matrix $A,
where op(A) = A, A^T, A^H for $TransA = $CblasNoTrans, $CblasTrans,
$CblasConjTrans (constant values coming from the CBLAS module). When $Uplo
is $CblasUpper then the upper triangle of $A is used, and when $Uplo is
$CblasLower then the lower triangle of $A is used. If $Diag is $CblasNonUnit
then the diagonal of the matrix is used, but if $Diag is $CblasUnit then the
diagonal elements of the matrix $A are taken as unity and are not
referenced. The function returns 0 if the operation suceeded, 1
otherwise.
- "gsl_blas_dtrsv($Uplo, $TransA, $Diag, $A, $x)" - This function
computes inv(op(A)) x for the vector $x, where op(A) = A, A^T, A^H for
$TransA = $CblasNoTrans, $CblasTrans, $CblasConjTrans (constant values
coming from the CBLAS module). When $Uplo is $CblasUpper then the upper
triangle of $A is used, and when $Uplo is $CblasLower then the lower
triangle of $A is used. If $Diag is $CblasNonUnit then the diagonal of the
matrix is used, but if $Diag is $CblasUnit then the diagonal elements of the
matrix $A are taken as unity and are not referenced. The function returns 0
if the operation suceeded, 1 otherwise.
- "gsl_blas_cgemv "
- "gsl_blas_ctrmv"
- "gsl_blas_ctrsv"
- "gsl_blas_zgemv "
- "gsl_blas_ztrmv"
- "gsl_blas_ztrsv"
- "gsl_blas_ssymv"
- "gsl_blas_sger "
- "gsl_blas_ssyr"
- "gsl_blas_ssyr2"
- "gsl_blas_dsymv"
- "gsl_blas_dger($alpha, $x, $y, $A)" - This function computes the
rank-1 update A = alpha x y^T + A of the matrix $A. $x and $y are vectors.
The function returns 0 if the operation suceeded, 1 otherwise.
- "gsl_blas_dsyr($Uplo, $alpha, $x, $A)" - This function computes
the symmetric rank-1 update A = \alpha x x^T + A of the symmetric matrix $A
and the vector $x. Since the matrix $A is symmetric only its upper half or
lower half need to be stored. When $Uplo is $CblasUpper then the upper
triangle and diagonal of $A are used, and when $Uplo is $CblasLower then the
lower triangle and diagonal of $A are used. The function returns 0 if the
operation suceeded, 1 otherwise.
- "gsl_blas_dsyr2($Uplo, $alpha, $x, $y, $A)" - This function
computes the symmetric rank-2 update A = \alpha x y^T + \alpha y x^T + A of
the symmetric matrix $A, the vector $x and vector $y. Since the matrix $A is
symmetric only its upper half or lower half need to be stored. When $Uplo is
$CblasUpper then the upper triangle and diagonal of $A are used, and when
$Uplo is $CblasLower then the lower triangle and diagonal of $A are
used.
- "gsl_blas_chemv"
- "gsl_blas_cgeru "
- "gsl_blas_cgerc"
- "gsl_blas_cher"
- "gsl_blas_cher2"
- "gsl_blas_zhemv "
- "gsl_blas_zgeru($alpha, $x, $y, $A)" - This function computes
the rank-1 update A = alpha x y^T + A of the complex matrix $A. $alpha is a
complex number and $x and $y are complex vectors. The function returns 0 if
the operation suceeded, 1 otherwise.
- "gsl_blas_zgerc"
- "gsl_blas_zher($Uplo, $alpha, $x, $A)" - This function computes
the hermitian rank-1 update A = \alpha x x^H + A of the hermitian matrix $A
and of the complex vector $x. Since the matrix $A is hermitian only its
upper half or lower half need to be stored. When $Uplo is $CblasUpper then
the upper triangle and diagonal of $A are used, and when $Uplo is
$CblasLower then the lower triangle and diagonal of $A are used. The
imaginary elements of the diagonal are automatically set to zero. The
function returns 0 if the operation suceeded, 1 otherwise.
- "gsl_blas_zher2 "
- "gsl_blas_sgemm"
- "gsl_blas_ssymm"
- "gsl_blas_ssyrk"
- "gsl_blas_ssyr2k "
- "gsl_blas_strmm"
- "gsl_blas_strsm"
- "gsl_blas_dgemm($TransA, $TransB, $alpha, $A, $B, $beta, $C)" -
This function computes the matrix-matrix product and sum C = \alpha op(A)
op(B) + \beta C where op(A) = A, A^T, A^H for $TransA = $CblasNoTrans,
$CblasTrans, $CblasConjTrans and similarly for the parameter $TransB. The
function returns 0 if the operation suceeded, 1 otherwise.
- "gsl_blas_dsymm($Side, $Uplo, $alpha, $A, $B, $beta, $C)" - This
function computes the matrix-matrix product and sum C = \alpha A B + \beta C
for $Side is $CblasLeft and C = \alpha B A + \beta C for $Side is
$CblasRight, where the matrix $A is symmetric. When $Uplo is $CblasUpper
then the upper triangle and diagonal of $A are used, and when $Uplo is
$CblasLower then the lower triangle and diagonal of $A are used. The
function returns 0 if the operation suceeded, 1 otherwise.
- "gsl_blas_dsyrk($Uplo, $Trans, $alpha, $A, $beta, $C)" - This
function computes a rank-k update of the symmetric matrix $C, C = \alpha A
A^T + \beta C when $Trans is $CblasNoTrans and C = \alpha A^T A + \beta C
when $Trans is $CblasTrans. Since the matrix $C is symmetric only its upper
half or lower half need to be stored. When $Uplo is $CblasUpper then the
upper triangle and diagonal of $C are used, and when $Uplo is $CblasLower
then the lower triangle and diagonal of $C are used. The function returns 0
if the operation suceeded, 1 otherwise.
- "gsl_blas_dsyr2k($Uplo, $Trans, $alpha, $A, $B, $beta, $C)" -
This function computes a rank-2k update of the symmetric matrix $C, C =
\alpha A B^T + \alpha B A^T + \beta C when $Trans is $CblasNoTrans and C =
\alpha A^T B + \alpha B^T A + \beta C when $Trans is $CblasTrans. Since the
matrix $C is symmetric only its upper half or lower half need to be stored.
When $Uplo is $CblasUpper then the upper triangle and diagonal of $C are
used, and when $Uplo is $CblasLower then the lower triangle and diagonal of
$C are used. The function returns 0 if the operation suceeded, 1
otherwise.
- "gsl_blas_dtrmm($Side, $Uplo, $TransA, $Diag, $alpha, $A, $B)" -
This function computes the matrix-matrix product B = \alpha op(A) B for
$Side is $CblasLeft and B = \alpha B op(A) for $Side is $CblasRight. The
matrix $A is triangular and op(A) = A, A^T, A^H for $TransA = $CblasNoTrans,
$CblasTrans, $CblasConjTrans. When $Uplo is $CblasUpper then the upper
triangle of $A is used, and when $Uplo is $CblasLower then the lower
triangle of $A is used. If $Diag is $CblasNonUnit then the diagonal of $A is
used, but if $Diag is $CblasUnit then the diagonal elements of the matrix $A
are taken as unity and are not referenced. The function returns 0 if the
operation suceeded, 1 otherwise.
- "gsl_blas_dtrsm($Side, $Uplo, $TransA, $Diag, $alpha, $A, $B)" -
This function computes the inverse-matrix matrix product B = \alpha
op(inv(A))B for $Side is $CblasLeft and B = \alpha B op(inv(A)) for $Side is
$CblasRight. The matrix $A is triangular and op(A) = A, A^T, A^H for $TransA
= $CblasNoTrans, $CblasTrans, $CblasConjTrans. When $Uplo is $CblasUpper
then the upper triangle of $A is used, and when $Uplo is $CblasLower then
the lower triangle of $A is used. If $Diag is $CblasNonUnit then the
diagonal of $A is used, but if $Diag is $CblasUnit then the diagonal
elements of the matrix $A are taken as unity and are not referenced. The
function returns 0 if the operation suceeded, 1 otherwise.
- "gsl_blas_cgemm"
- "gsl_blas_csymm"
- "gsl_blas_csyrk"
- "gsl_blas_csyr2k "
- "gsl_blas_ctrmm"
- "gsl_blas_ctrsm"
- "gsl_blas_zgemm($TransA, $TransB, $alpha, $A, $B, $beta, $C)" -
This function computes the matrix-matrix product and sum C = \alpha op(A)
op(B) + \beta C where op(A) = A, A^T, A^H for $TransA = $CblasNoTrans,
$CblasTrans, $CblasConjTrans and similarly for the parameter $TransB. The
function returns 0 if the operation suceeded, 1 otherwise. $A, $B and $C are
complex matrices
- "gsl_blas_zsymm($Side, $Uplo, $alpha, $A, $B, $beta, $C)" - This
function computes the matrix-matrix product and sum C = \alpha A B + \beta C
for $Side is $CblasLeft and C = \alpha B A + \beta C for $Side is
$CblasRight, where the matrix $A is symmetric. When $Uplo is $CblasUpper
then the upper triangle and diagonal of $A are used, and when $Uplo is
$CblasLower then the lower triangle and diagonal of $A are used. $A, $B and
$C are complex matrices. The function returns 0 if the operation suceeded, 1
otherwise.
- "gsl_blas_zsyrk($Uplo, $Trans, $alpha, $A, $beta, $C)" - This
function computes a rank-k update of the symmetric complex matrix $C, C =
\alpha A A^T + \beta C when $Trans is $CblasNoTrans and C = \alpha A^T A +
\beta C when $Trans is $CblasTrans. Since the matrix $C is symmetric only
its upper half or lower half need to be stored. When $Uplo is $CblasUpper
then the upper triangle and diagonal of $C are used, and when $Uplo is
$CblasLower then the lower triangle and diagonal of $C are used. The
function returns 0 if the operation suceeded, 1 otherwise.
- "gsl_blas_zsyr2k($Uplo, $Trans, $alpha, $A, $B, $beta, $C)" -
This function computes a rank-2k update of the symmetric matrix $C, C =
\alpha A B^T + \alpha B A^T + \beta C when $Trans is $CblasNoTrans and C =
\alpha A^T B + \alpha B^T A + \beta C when $Trans is $CblasTrans. Since the
matrix $C is symmetric only its upper half or lower half need to be stored.
When $Uplo is $CblasUpper then the upper triangle and diagonal of $C are
used, and when $Uplo is $CblasLower then the lower triangle and diagonal of
$C are used. The function returns 0 if the operation suceeded, 1 otherwise.
$A, $B and $C are complex matrices and $beta is a complex number.
- "gsl_blas_ztrmm($Side, $Uplo, $TransA, $Diag, $alpha, $A, $B)" -
This function computes the matrix-matrix product B = \alpha op(A) B for
$Side is $CblasLeft and B = \alpha B op(A) for $Side is $CblasRight. The
matrix $A is triangular and op(A) = A, A^T, A^H for $TransA = $CblasNoTrans,
$CblasTrans, $CblasConjTrans. When $Uplo is $CblasUpper then the upper
triangle of $A is used, and when $Uplo is $CblasLower then the lower
triangle of $A is used. If $Diag is $CblasNonUnit then the diagonal of $A is
used, but if $Diag is $CblasUnit then the diagonal elements of the matrix $A
are taken as unity and are not referenced. The function returns 0 if the
operation suceeded, 1 otherwise. $A and $B are complex matrices and $alpha
is a complex number.
- "gsl_blas_ztrsm($Side, $Uplo, $TransA, $Diag, $alpha, $A, $B)" -
This function computes the inverse-matrix matrix product B = \alpha
op(inv(A))B for $Side is $CblasLeft and B = \alpha B op(inv(A)) for $Side is
$CblasRight. The matrix $A is triangular and op(A) = A, A^T, A^H for $TransA
= $CblasNoTrans, $CblasTrans, $CblasConjTrans. When $Uplo is $CblasUpper
then the upper triangle of $A is used, and when $Uplo is $CblasLower then
the lower triangle of $A is used. If $Diag is $CblasNonUnit then the
diagonal of $A is used, but if $Diag is $CblasUnit then the diagonal
elements of the matrix $A are taken as unity and are not referenced. The
function returns 0 if the operation suceeded, 1 otherwise. $A and $B are
complex matrices and $alpha is a complex number.
- "gsl_blas_chemm"
- "gsl_blas_cherk"
- "gsl_blas_cher2k"
- "gsl_blas_zhemm($Side, $Uplo, $alpha, $A, $B, $beta, $C)" - This
function computes the matrix-matrix product and sum C = \alpha A B + \beta C
for $Side is $CblasLeft and C = \alpha B A + \beta C for $Side is
$CblasRight, where the matrix $A is hermitian. When Uplo is CblasUpper then
the upper triangle and diagonal of A are used, and when Uplo is CblasLower
then the lower triangle and diagonal of A are used. The imaginary elements
of the diagonal are automatically set to zero.
- "gsl_blas_zherk($Uplo, $Trans, $alpha, $A, $beta, $C)" - This
function computes a rank-k update of the hermitian matrix $C, C = \alpha A
A^H + \beta C when $Trans is $CblasNoTrans and C = \alpha A^H A + \beta C
when $Trans is $CblasTrans. Since the matrix $C is hermitian only its upper
half or lower half need to be stored. When $Uplo is $CblasUpper then the
upper triangle and diagonal of $C are used, and when $Uplo is $CblasLower
then the lower triangle and diagonal of $C are used. The imaginary elements
of the diagonal are automatically set to zero. The function returns 0 if the
operation suceeded, 1 otherwise. $A, $B and $C are complex matrices and
$alpha and $beta are complex numbers.
- "gsl_blas_zher2k($Uplo, $Trans, $alpha, $A, $B, $beta, $C)" -
This function computes a rank-2k update of the hermitian matrix $C, C =
\alpha A B^H + \alpha^* B A^H + \beta C when $Trans is $CblasNoTrans and C =
\alpha A^H B + \alpha^* B^H A + \beta C when $Trans is $CblasConjTrans.
Since the matrix $C is hermitian only its upper half or lower half need to
be stored. When $Uplo is $CblasUpper then the upper triangle and diagonal of
$C are used, and when $Uplo is $CblasLower then the lower triangle and
diagonal of $C are used. The imaginary elements of the diagonal are
automatically set to zero. The function returns 0 if the operation suceeded,
1 otherwise.
You have to add the functions you want to use inside the qw
/put_function_here /. You can also write use Math::GSL::BLAS qw/:all/ to use
all available functions of the module. Other tags are also available, here
is a complete list of all tags for this module :
- "level1"
- "level2"
- "level3"
For more informations on the functions, we refer you to the GSL
official documentation:
<http://www.gnu.org/software/gsl/manual/html_node/>
Jonathan "Duke" Leto <jonathan@leto.net> and Thierry Moisan
<thierry.moisan@gmail.com>
Copyright (C) 2008-2021 Jonathan "Duke" Leto and Thierry Moisan
This program is free software; you can redistribute it and/or
modify it under the same terms as Perl itself.
Visit the GSP FreeBSD Man Page Interface. Output converted with ManDoc. |