File: BLAS.pod

package info (click to toggle)
libmath-gsl-perl 0.45-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 192,156 kB
  • sloc: ansic: 895,524; perl: 24,682; makefile: 12
file content (383 lines) | stat: -rw-r--r-- 21,334 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
%perlcode %{

@EXPORT_OK_level1 = qw/
                        gsl_blas_sdsdot gsl_blas_dsdot gsl_blas_sdot gsl_blas_ddot
                        gsl_blas_cdotu gsl_blas_cdotc gsl_blas_zdotu gsl_blas_zdotc
                        gsl_blas_snrm2 gsl_blas_sasum gsl_blas_dnrm2 gsl_blas_dasum
                        gsl_blas_scnrm2 gsl_blas_scasum gsl_blas_dznrm2 gsl_blas_dzasum
                        gsl_blas_isamax gsl_blas_idamax gsl_blas_icamax gsl_blas_izamax
                        gsl_blas_sswap gsl_blas_scopy gsl_blas_saxpy gsl_blas_dswap
                        gsl_blas_dcopy gsl_blas_daxpy 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 gsl_blas_drotm
                        gsl_blas_sscal gsl_blas_dscal gsl_blas_cscal gsl_blas_zscal
                        gsl_blas_csscal gsl_blas_zdscal
                    /;
@EXPORT_OK_level2 = qw/
                        gsl_blas_sgemv gsl_blas_strmv
                        gsl_blas_strsv gsl_blas_dgemv gsl_blas_dtrmv gsl_blas_dtrsv
                        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
                        gsl_blas_dsyr gsl_blas_dsyr2 gsl_blas_chemv gsl_blas_cgeru
                        gsl_blas_cgerc gsl_blas_cher gsl_blas_cher2 gsl_blas_zhemv
                        gsl_blas_zgeru gsl_blas_zgerc gsl_blas_zher gsl_blas_zher2
                    /;

@EXPORT_OK_level3 = qw/
                        gsl_blas_sgemm gsl_blas_ssymm gsl_blas_ssyrk gsl_blas_ssyr2k
                        gsl_blas_strmm gsl_blas_strsm gsl_blas_dgemm gsl_blas_dsymm
                        gsl_blas_dsyrk gsl_blas_dsyr2k gsl_blas_dtrmm gsl_blas_dtrsm
                        gsl_blas_cgemm gsl_blas_csymm gsl_blas_csyrk gsl_blas_csyr2k
                        gsl_blas_ctrmm gsl_blas_ctrsm gsl_blas_zgemm gsl_blas_zsymm
                        gsl_blas_zsyrk gsl_blas_zsyr2k gsl_blas_ztrmm gsl_blas_ztrsm
                        gsl_blas_chemm gsl_blas_cherk gsl_blas_cher2k gsl_blas_zhemm
                        gsl_blas_zherk gsl_blas_zher2k
                    /;
@EXPORT_OK = (@EXPORT_OK_level1, @EXPORT_OK_level2, @EXPORT_OK_level3);
%EXPORT_TAGS = (
                all    => [ @EXPORT_OK ],
                level1 => [ @EXPORT_OK_level1 ],
                level2 => [ @EXPORT_OK_level2 ],
                level3 => [ @EXPORT_OK_level3 ],
               );

__END__

=encoding utf8

=head1 NAME

Math::GSL::BLAS - Basic Linear Algebra Subprograms

=head1 SYNOPSIS

    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";

=head1 DESCRIPTION

The functions of this module are divised into 3 levels:

=head2 Level 1 - Vector operations

=over 3

=item C<gsl_blas_sdsdot>

=item C<gsl_blas_dsdot>

=item C<gsl_blas_sdot>

=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 succeeded, 1
otherwise and the second value is the result of the computation.

=item C<gsl_blas_cdotu>

=item C<gsl_blas_cdotc>

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise.

=item C<gsl_blas_snrm2>
=item C<gsl_blas_sasum>

=item C<gsl_blas_dnrm2($x)>

This function computes the Euclidean norm

    ||x||_2 = \sqrt {\sum x_i^2}

of the vector $x.

=item C<gsl_blas_dasum($x)>

This function computes the absolute sum \sum |x_i| of the elements of the vector $x.

=item C<gsl_blas_scnrm2>

=item C<gsl_blas_scasum>

=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)}.

=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)|.

=item C<gsl_blas_isamax>

=item C<gsl_blas_idamax>

=item C<gsl_blas_icamax>

=item C<gsl_blas_izamax >

=item C<gsl_blas_sswap>

=item C<gsl_blas_scopy>

=item C<gsl_blas_saxpy>

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise.

=item C<gsl_blas_daxpy($alpha, $x, $y)>

These functions compute the sum $y = $alpha * $x + $y for the vectors $x and $y.

=item C<gsl_blas_cswap>

=item C<gsl_blas_ccopy >

=item C<gsl_blas_caxpy>

=item C<gsl_blas_zswap>

=item C<gsl_blas_zcopy>

=item C<gsl_blas_zaxpy >

=item C<gsl_blas_srotg>

=item C<gsl_blas_srotmg>

=item C<gsl_blas_srot>

=item C<gsl_blas_srotm >

=item C<gsl_blas_drotg>

=item C<gsl_blas_drotmg>

=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.

=item C<gsl_blas_drotm >

=item C<gsl_blas_sscal>

=item C<gsl_blas_dscal($alpha, $x)>

This function rescales the vector $x by the multiplicative factor $alpha.

=item C<gsl_blas_cscal>

=item C<gsl_blas_zscal >

=item C<gsl_blas_csscal>

=item C<gsl_blas_zdscal>

=back

=head2 Level 2 - Matrix-vector operations

=over 3

=item C<gsl_blas_sgemv>

=item C<gsl_blas_strmv >

=item C<gsl_blas_strsv>

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise.

=item C<gsl_blas_cgemv >

=item C<gsl_blas_ctrmv>

=item C<gsl_blas_ctrsv>

=item C<gsl_blas_zgemv >

=item C<gsl_blas_ztrmv>

=item C<gsl_blas_ztrsv>

=item C<gsl_blas_ssymv>

=item C<gsl_blas_sger >

=item C<gsl_blas_ssyr>

=item C<gsl_blas_ssyr2>

=item C<gsl_blas_dsymv>

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise.

=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.

=item C<gsl_blas_chemv>

=item C<gsl_blas_cgeru >

=item C<gsl_blas_cgerc>

=item C<gsl_blas_cher>

=item C<gsl_blas_cher2>

=item C<gsl_blas_zhemv >

=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 succeeded, 1 otherwise.

=item C<gsl_blas_zgerc>

=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 succeeded, 1 otherwise.


=item C<gsl_blas_zher2 >

=back

=head2 Level 3 - Matrix-matrix operations

=over 3

=item C<gsl_blas_sgemm>

=item C<gsl_blas_ssymm>

=item C<gsl_blas_ssyrk>

=item C<gsl_blas_ssyr2k >

=item C<gsl_blas_strmm>

=item C<gsl_blas_strsm>

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise.

=item C<gsl_blas_cgemm>

=item C<gsl_blas_csymm>

=item C<gsl_blas_csyrk>

=item C<gsl_blas_csyr2k >

=item C<gsl_blas_ctrmm>

=item C<gsl_blas_ctrsm>

=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 succeeded, 1 otherwise. $A, $B and $C are complex matrices

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise.

=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 succeeded, 1 otherwise. $A, $B and $C are complex matrices and $beta is a complex number.

=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 succeeded, 1 otherwise. $A and $B are complex matrices and $alpha is a complex number.

=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 succeeded, 1 otherwise. $A and $B are complex matrices and $alpha is a complex number.

=item C<gsl_blas_chemm>

=item C<gsl_blas_cherk>

=item C<gsl_blas_cher2k>

=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.

=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 succeeded, 1 otherwise. $A, $B and $C are complex matrices and $alpha and $beta are complex numbers.

=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 succeeded, 1 otherwise.

=back

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 :

=over 3

=item C<level1>

=item C<level2>

=item C<level3>

=back

For more information on the functions, we refer you to the GSL official documentation: L<http://www.gnu.org/software/gsl/manual/html_node/>

=head1 AUTHORS

Jonathan "Duke" Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>

=head1 COPYRIGHT AND LICENSE

Copyright (C) 2008-2024 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.

=cut

%}