Added myself to the CREDITS file (after Jonathan Leto reminded me about it)
[Math-GSL.git] / pod / BLAS.pod
blobd1ed7cb855aedebe1483531c34f7302803318c86
1 %perlcode %{
3 @EXPORT_OK_level1 = qw/
4                         gsl_blas_sdsdot gsl_blas_dsdot gsl_blas_sdot gsl_blas_ddot 
5                         gsl_blas_cdotu gsl_blas_cdotc gsl_blas_zdotu gsl_blas_zdotc 
6                         gsl_blas_snrm2 gsl_blas_sasum gsl_blas_dnrm2 gsl_blas_dasum 
7                         gsl_blas_scnrm2 gsl_blas_scasum gsl_blas_dznrm2 gsl_blas_dzasum 
8                         gsl_blas_isamax gsl_blas_idamax gsl_blas_icamax gsl_blas_izamax 
9                         gsl_blas_sswap gsl_blas_scopy gsl_blas_saxpy gsl_blas_dswap 
10                         gsl_blas_dcopy gsl_blas_daxpy gsl_blas_cswap gsl_blas_ccopy 
11                         gsl_blas_caxpy gsl_blas_zswap gsl_blas_zcopy gsl_blas_zaxpy 
12                         gsl_blas_srotg gsl_blas_srotmg gsl_blas_srot gsl_blas_srotm 
13                         gsl_blas_drotg gsl_blas_drotmg gsl_blas_drot gsl_blas_drotm 
14                         gsl_blas_sscal gsl_blas_dscal gsl_blas_cscal gsl_blas_zscal 
15                         gsl_blas_csscal gsl_blas_zdscal
16                     /;
17 @EXPORT_OK_level2 = qw/
18                         gsl_blas_sgemv gsl_blas_strmv 
19                         gsl_blas_strsv gsl_blas_dgemv gsl_blas_dtrmv gsl_blas_dtrsv 
20                         gsl_blas_cgemv gsl_blas_ctrmv gsl_blas_ctrsv gsl_blas_zgemv 
21                         gsl_blas_ztrmv gsl_blas_ztrsv gsl_blas_ssymv gsl_blas_sger 
22                         gsl_blas_ssyr gsl_blas_ssyr2 gsl_blas_dsymv gsl_blas_dger 
23                         gsl_blas_dsyr gsl_blas_dsyr2 gsl_blas_chemv gsl_blas_cgeru 
24                         gsl_blas_cgerc gsl_blas_cher gsl_blas_cher2 gsl_blas_zhemv 
25                         gsl_blas_zgeru gsl_blas_zgerc gsl_blas_zher gsl_blas_zher2 
26                     /;
28 @EXPORT_OK_level3 = qw/
29                         gsl_blas_sgemm gsl_blas_ssymm gsl_blas_ssyrk gsl_blas_ssyr2k 
30                         gsl_blas_strmm gsl_blas_strsm gsl_blas_dgemm gsl_blas_dsymm 
31                         gsl_blas_dsyrk gsl_blas_dsyr2k gsl_blas_dtrmm gsl_blas_dtrsm 
32                         gsl_blas_cgemm gsl_blas_csymm gsl_blas_csyrk gsl_blas_csyr2k 
33                         gsl_blas_ctrmm gsl_blas_ctrsm gsl_blas_zgemm gsl_blas_zsymm 
34                         gsl_blas_zsyrk gsl_blas_zsyr2k gsl_blas_ztrmm gsl_blas_ztrsm 
35                         gsl_blas_chemm gsl_blas_cherk gsl_blas_cher2k gsl_blas_zhemm 
36                         gsl_blas_zherk gsl_blas_zher2k 
37                     /;
38 @EXPORT_OK = (@EXPORT_OK_level1, @EXPORT_OK_level2, @EXPORT_OK_level3);
39 %EXPORT_TAGS = (
40                 all    => [ @EXPORT_OK ],
41                 level1 => [ @EXPORT_OK_level1 ],  
42                 level2 => [ @EXPORT_OK_level2 ],  
43                 level3 => [ @EXPORT_OK_level3 ],  
44                );
46 __END__
48 =head1 NAME
50 Math::GSL::BLAS - Basic Linear Algebra Subprograms
52 =head1 SYNOPSIS
54 use Math::GSL::BLAS qw/:all/;
56 =head1 DESCRIPTION
58 The functions of this module are divised into 3 levels:
60 =head2 Level 1 - Vector operations
62 =over 3  
64 =item C<gsl_blas_sdsdot>
66 =item C<gsl_blas_dsdot>
68 =item C<gsl_blas_sdot>
70 =item C<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. 
72 =item C<gsl_blas_cdotu>
74 =item C<gsl_blas_cdotc>
76 =item C<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.
78 =item C<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.
80 =item C<gsl_blas_snrm2> 
81 =item C<gsl_blas_sasum>
83 =item C<gsl_blas_dnrm2($x)> - This function computes the Euclidean norm ||x||_2 = \sqrt {\sum x_i^2} of the vector $x. 
85 =item C<gsl_blas_dasum($x)> - This function computes the absolute sum \sum |x_i| of the elements of the vector $x. 
87 =item C<gsl_blas_scnrm2>
89 =item C<gsl_blas_scasum>
91 =item C<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)}.
93 =item C<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)|. 
95 =item C<gsl_blas_isamax>
97 =item C<gsl_blas_idamax>
99 =item C<gsl_blas_icamax>
101 =item C<gsl_blas_izamax >
103 =item C<gsl_blas_sswap>
105 =item C<gsl_blas_scopy>
107 =item C<gsl_blas_saxpy>
109 =item C<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.
111 =item C<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.
113 =item C<gsl_blas_daxpy($alpha, $x, $y)> - These functions compute the sum $y = $alpha * $x + $y for the vectors $x and $y. 
115 =item C<gsl_blas_cswap>
117 =item C<gsl_blas_ccopy >
119 =item C<gsl_blas_caxpy>
121 =item C<gsl_blas_zswap>
123 =item C<gsl_blas_zcopy>
125 =item C<gsl_blas_zaxpy >
127 =item C<gsl_blas_srotg>
129 =item C<gsl_blas_srotmg>
131 =item C<gsl_blas_srot>
133 =item C<gsl_blas_srotm >
135 =item C<gsl_blas_drotg>
137 =item C<gsl_blas_drotmg>
139 =item C<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. 
141 =item C<gsl_blas_drotm >
143 =item C<gsl_blas_sscal>
145 =item C<gsl_blas_dscal($alpha, $x)> - This function rescales the vector $x by the multiplicative factor $alpha.
147 =item C<gsl_blas_cscal>
149 =item C<gsl_blas_zscal >
151 =item C<gsl_blas_csscal>
153 =item C<gsl_blas_zdscal>
155 =back
157 =head2 Level 2 - Matrix-vector operations
159 =over 3 
161 =item C<gsl_blas_sgemv>
163 =item C<gsl_blas_strmv >
165 =item C<gsl_blas_strsv>
167 =item C<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.
169 =item C<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.
171 =item C<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.
173 =item C<gsl_blas_cgemv >
175 =item C<gsl_blas_ctrmv>
177 =item C<gsl_blas_ctrsv>
179 =item C<gsl_blas_zgemv >
181 =item C<gsl_blas_ztrmv>
183 =item C<gsl_blas_ztrsv>
185 =item C<gsl_blas_ssymv>
187 =item C<gsl_blas_sger >
189 =item C<gsl_blas_ssyr>
191 =item C<gsl_blas_ssyr2>
193 =item C<gsl_blas_dsymv>
195 =item C<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.
197 =item C<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.
199 =item C<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.
201 =item C<gsl_blas_chemv>
203 =item C<gsl_blas_cgeru >
205 =item C<gsl_blas_cgerc>
207 =item C<gsl_blas_cher>
209 =item C<gsl_blas_cher2>
211 =item C<gsl_blas_zhemv >
213 =item C<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.
215 =item C<gsl_blas_zgerc>
217 =item C<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.
220 =item C<gsl_blas_zher2 >
222 =back
224 =head2 Level 3 - Matrix-matrix operations
226 =over 3 
228 =item C<gsl_blas_sgemm>
230 =item C<gsl_blas_ssymm>
232 =item C<gsl_blas_ssyrk>
234 =item C<gsl_blas_ssyr2k >
236 =item C<gsl_blas_strmm>
238 =item C<gsl_blas_strsm>
240 =item C<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.
242 =item C<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.
244 =item C<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.
246 =item C<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.
248 =item C<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.
250 =item C<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.
252 =item C<gsl_blas_cgemm>
254 =item C<gsl_blas_csymm>
256 =item C<gsl_blas_csyrk>
258 =item C<gsl_blas_csyr2k >
260 =item C<gsl_blas_ctrmm>
262 =item C<gsl_blas_ctrsm>
264 =item C<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
266 =item C<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.
268 =item C<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.
270 =item C<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.
272 =item C<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.
274 =item C<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.
276 =item C<gsl_blas_chemm>
278 =item C<gsl_blas_cherk>
280 =item C<gsl_blas_cher2k>
282 =item C<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. 
284 =item C<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.
286 =item C<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.
288 =back
290 You have to add the functions you want to use inside the qw /put_funtion_here /. 
291 You can also write use Math::GSL::BLAS qw/:all/ to use all avaible functions of the module. 
292 Other tags are also avaible, here is a complete list of all tags for this module :
294 =over 3
296 =item C<level1>
298 =item C<level2>
300 =item C<level3> 
302 =back
304 For more informations on the functions, we refer you to the GSL offcial documentation: L<http://www.gnu.org/software/gsl/manual/html_node/>
306 Tip : search on google: site:http://www.gnu.org/software/gsl/manual/html_node/ name_of_the_function_you_want
309 =head1 EXAMPLES
311  This example shows how to do a matrix-matrix product of double numbers :
313  use Math::GSL::Matrix qw/:all/;
314  use Math::GSL::BLAS qw/:all/;
315  my $A = Math::GSL::Matrix->new(2,2);
316  $A->set_row(0, [1, 4]);
317    ->set_row(1, [3, 2]);
318  my $B = Math::GSL::Matrix->new(2,2);
319  $B->set_row(0, [2, 1]);
320    ->set_row(1, [5,3]);
321  my $C = Math::GSL::Matrix->new(2,2);
322  gsl_matrix_set_zero($C->raw);
323  gsl_blas_dgemm($CblasNoTrans, $CblasNoTrans, 1, $A->raw, $B->raw, 1, $C->raw);
324  my @got = $C->row(0)->as_list;
325  print "The resulting matrix is: \n[";
326  print "$got[0]  $got[1]\n";
327  @got = $C->row(1)->as_list;
328  print "$got[0]  $got[1] ]\n";
331  This example shows how to compute the scalar product of two vectors :
333  use Math::GSL::Vector qw/:all/;
334  use Math::GSL::CBLAS qw/:all/;
335  use Math::GSL::BLAS qw/:all/;
336  my $vec1 = Math::GSL::Vector->new([1,2,3,4,5]);
337  my $vec2 = Math::GSL::Vector->new([5,4,3,2,1]);
338  my ($status, $result) = gsl_blas_ddot($vec1->raw, $vec2->raw);
339  if($status == 0) { 
340  print "The function has succeeded. \n";
342  print "The result of the vector multiplication is $result. \n";
344 =head1 AUTHORS
346 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
348 =head1 COPYRIGHT AND LICENSE
350 Copyright (C) 2008-2009 Jonathan Leto and Thierry Moisan
352 This program is free software; you can redistribute it and/or modify it
353 under the same terms as Perl itself.
355 =cut