Refactor inclusion of system.i to one place
[Math-GSL.git] / t / BLAS.t
bloba8afe30babf8b3a6432df618ae9a97c34c7d66ac
1 package Math::GSL::BLAS::Test;
2 use base q{Test::Class};
3 use Test::More tests => 99;
4 use Math::GSL          qw/:all/;
5 use Math::GSL::BLAS    qw/:all/;
6 use Math::GSL::Vector  qw/:all/;
7 use Math::GSL::Complex qw/:all/;
8 use Math::GSL::Matrix  qw/:all/;
9 use Math::GSL::CBLAS   qw/:all/;
10 use Math::GSL::Test    qw/:all/;
11 use Math::GSL::Errno   qw/:all/;
12 use Math::GSL::MatrixComplex  qw/:all/;
13 use Math::GSL::VectorComplex  qw/:all/;
14 use Data::Dumper;
15 use strict;
17 BEGIN { gsl_set_error_handler_off(); }
19 sub make_fixture : Test(setup) {
22 sub teardown : Test(teardown) {
25 sub GSL_BLAS_DDOT : Tests {
26     my $vec1 = Math::GSL::Vector->new([1,2,3,4,5]);
27     my $vec2 = Math::GSL::Vector->new([5,4,3,2,1]);
28     my ($x, $result) = gsl_blas_ddot($vec1->raw, $vec2->raw);
29     ok_status($x);
30     is_similar($result,35);
33 sub GSL_BLAS_ZDOTU : Tests {
34     my $vec1 = gsl_vector_complex_alloc(2);
35     my $vec2 = gsl_vector_complex_alloc(2);
36     my $c = gsl_complex_rect(2,1);
37     gsl_vector_complex_set($vec1,0,$c); 
38     gsl_vector_complex_set($vec2,0,$c); 
39     $c = gsl_complex_rect(1,1);
40     gsl_vector_complex_set($vec1,1,$c); 
41     gsl_vector_complex_set($vec2,1,$c);
42     ok_status(gsl_blas_zdotu($vec1, $vec2, $c));
43     is_similar([ gsl_parts($c) ], [3,6]); 
46 sub GSL_BLAS_ZDOTC : Tests {
47     my $vec1 = gsl_vector_complex_alloc(2);
48     my $vec2 = gsl_vector_complex_alloc(2);
49     my $c = gsl_complex_rect(2,1);
50     gsl_vector_complex_set($vec1,0,$c); 
51     gsl_vector_complex_set($vec2,0,$c); 
52     $c = gsl_complex_rect(1,1);
53     gsl_vector_complex_set($vec1,1,$c); 
54     gsl_vector_complex_set($vec2,1,$c);
55     ok_status(gsl_blas_zdotc($vec1, $vec2, $c));
56     is_similar([gsl_parts($c)], [7,0]); 
59 sub GSL_BLAS_DNRM2 : Tests {
60     my $vec = Math::GSL::Vector->new([3,4]);
61     is_similar(gsl_blas_dnrm2($vec->raw), 5);
63   
65 sub GSL_BLAS_DZNRM2 : Tests {
66     my $vec = gsl_vector_complex_alloc(2);
67     my $c = gsl_complex_rect(2,1);
68     gsl_vector_complex_set($vec,0,$c); 
69     $c = gsl_complex_rect(1,1);
70     gsl_vector_complex_set($vec,1,$c); 
71     is_similar([gsl_blas_dznrm2($vec)], [sqrt(7)]);
74 sub GSL_BLAS_DASUM : Tests {
75     my $vec = Math::GSL::Vector->new([2,-3,4]);
76     is_similar(gsl_blas_dasum($vec->raw), 9);
79 sub GSL_BLAS_DZASUM : Tests {
80     my $vec = gsl_vector_complex_alloc(2);
81     my $c = gsl_complex_rect(2,1);
82     gsl_vector_complex_set($vec,0,$c); 
83     $c = gsl_complex_rect(1,1);
84     gsl_vector_complex_set($vec,1,$c); 
85     is_similar(gsl_blas_dzasum($vec), 5);
88 sub GSL_BLAS_DSWAP : Tests {
89     my $vec1 = Math::GSL::Vector->new([0,1,2]);
90     my $vec2 = Math::GSL::Vector->new([2,1,0]);
91     gsl_blas_dswap($vec1->raw, $vec2->raw);
92     ok_similar( [0 .. 2], [ $vec2->as_list ] );
93     ok_similar( [2, 1,0], [ $vec1->as_list ] );
96 sub GSL_BLAS_ZSWAP : Tests { 
97     my $vec1 = gsl_vector_complex_alloc(2);
98     my $vec2 = gsl_vector_complex_alloc(2);
99     my $c = gsl_complex_rect(5,4);
100     gsl_vector_complex_set($vec1,0,$c); 
101     $c = gsl_complex_rect(2,2);
102     gsl_vector_complex_set($vec1,1, $c);
103     $c = gsl_complex_rect(3,3);
104     gsl_vector_complex_set($vec2,0, $c);
105     $c = gsl_complex_rect(1,1); 
106     gsl_vector_complex_set($vec2,1, $c);
108     ok_status(gsl_blas_zswap($vec1, $vec2));
109     $c = gsl_vector_complex_get($vec1,0);
110     ok( defined $c,"gsl_vector_complex_get");  
112     # goes boom
113     #is_similar( [gsl_parts($c) ], [ 3,3 ] );
116 sub GSL_BLAS_DCOPY : Tests {
117     my $vec1 = Math::GSL::Vector->new([0,1,2]);
118     my $vec2 = Math::GSL::Vector->new(3);
119     ok_status(gsl_blas_dcopy($vec1->raw, $vec2->raw));
120     ok_similar( [ $vec2->as_list ], [ 0 .. 2] );
123 sub GSL_BLAS_DAXPY : Tests { 
124     my $vec1 = Math::GSL::Vector->new([0,1,2]);
125     my $vec2 = Math::GSL::Vector->new([2,3,4]);
126     ok_status(gsl_blas_daxpy(2,$vec1->raw, $vec2->raw));
127     is_similar( [ $vec2->as_list ], [ 2, 5, 8 ] );
130 sub GSL_BLAS_DSCAL : Tests {
131     my $vec = Math::GSL::Vector->new([0,1,2]);
132     gsl_blas_dscal(4, $vec->raw);
133     is_similar( [ $vec->as_list ], [0,4,8] );
136 sub GSL_BLAS_DROT : Tests {
137     my $x = Math::GSL::Vector->new([1,2,3]);
138     my $y = Math::GSL::Vector->new([0,1,2]);
139     ok_status(gsl_blas_drot($x->raw,$y->raw,2,3));
140     ok_similar( [$x->as_list], [ 2,7,12], 'first vector');
141     ok_similar( [$y->as_list], [-3,-4,-5], 'second vector');
144 sub GSL_BLAS_DGER : Tests { 
145  my $x = Math::GSL::Vector->new([1,2,3]);
146  my $y = Math::GSL::Vector->new([0,1,2]);
147  my $A = Math::GSL::Matrix->new(3,3); 
148  gsl_matrix_set_zero($A->raw);
149  ok_status(gsl_blas_dger(2, $x->raw, $y->raw, $A->raw));
150  ok_similar([$A->row(0)->as_list], [0,2,4]);
151  ok_similar([$A->row(1)->as_list], [0,4,8]);
152  ok_similar([$A->row(2)->as_list], [0,6,12]);
155 sub GSL_BLAS_ZGERU : Tests {
156     my $x = gsl_vector_complex_alloc(2);
157     my $y = gsl_vector_complex_alloc(2);
158     my $A = gsl_matrix_complex_alloc(2,2);
159     my $alpha = gsl_complex_rect(2,2);
160     gsl_vector_complex_set($x, 0, $alpha);
161     $alpha = gsl_complex_rect(1,2);
162     gsl_vector_complex_set($x, 1, $alpha);
163     gsl_vector_complex_set($y, 0, $alpha);
164     $alpha = gsl_complex_rect(3,2);
165     gsl_vector_complex_set($y, 1, $alpha);
166     $alpha = gsl_complex_rect(0,0);
167     for (my $line=0; $line<2; $line++) {
168     map { gsl_matrix_complex_set($A, $line, $_, $alpha) } (0..1); }
169     $alpha = gsl_complex_rect(1,0);
170     ok_status(gsl_blas_zgeru($alpha, $x, $y, $A));
171     
172     $alpha= gsl_matrix_complex_get($A, 0,0);
173     ok_similar([gsl_parts($alpha)], [-2, 6]);
174     $alpha= gsl_matrix_complex_get($A, 1,0);
175     ok_similar([gsl_parts($alpha)], [-3, 4]);
176     $alpha= gsl_matrix_complex_get($A, 1,0);
177     ok_similar([gsl_parts($alpha)], [-3, 4]);
178     $alpha= gsl_matrix_complex_get($A, 0,1); 
179     ok_similar([gsl_parts($alpha)], [2, 10]);
180     $alpha= gsl_matrix_complex_get($A, 1,1); 
181     ok_similar([gsl_parts($alpha)], [-1, 8]);
184 sub GSL_BLAS_DGEMV : Tests {
185     my $x = Math::GSL::Vector->new([1,2,3]);
186     my $y = Math::GSL::Vector->new([0,1,2]);
187     my $A = Math::GSL::Matrix->new(3,3);
188     gsl_matrix_set_identity($A->raw);
189     ok_status(gsl_blas_dgemv($CblasNoTrans, 2, $A->raw, $x->raw,2, $y->raw));
190     ok_similar( [ $y->as_list ], [  2, 6, 10 ] );
193 sub GSL_BLAS_DTRMV : Tests {
194     my $x = Math::GSL::Vector->new([1,2,3]);
195     my $A = Math::GSL::Matrix->new(3,3);
196     gsl_matrix_set($A->raw, 0,0,3);
197     gsl_matrix_set($A->raw, 1,1,3);
198     gsl_matrix_set($A->raw, 2,2,3);
199     gsl_matrix_set($A->raw, 0,1,2);
200     gsl_matrix_set($A->raw, 0,2,3);
201     gsl_matrix_set($A->raw, 1,2,4);
202     ok_status(gsl_blas_dtrmv($CblasLower, $CblasNoTrans, $CblasNonUnit, $A->raw, $x->raw));
203     is_similar( [ $x->as_list ], [ 3, 6, 9 ] );
206 sub GSL_BLAS_DTRSV : Tests {
207     my $x = Math::GSL::Vector->new([40,40,40,40]);
208     my $A = Math::GSL::Matrix->new(4,4);
209     map { gsl_matrix_set($A->raw, $_,0,$_+1); } (0..3);
210     map { gsl_matrix_set($A->raw, $_,1,$_+2); } (0..2);
211     gsl_matrix_set($A->raw, 3,1,1);
212     map { gsl_matrix_set($A->raw, $_,2,$_+3); } (0..1);
213     map { gsl_matrix_set($A->raw, $_,2,$_-1); } (2..3);
214     map { gsl_matrix_set($A->raw, $_,3,$_); } (1..3);
215     gsl_matrix_set($A->raw, 0,3,4);
216     ok_status(gsl_blas_dtrsv($CblasLower, $CblasNoTrans, $CblasNonUnit, $A->raw, $x->raw));
217     ok_similar([$x->as_list], [40,-40/3,-80/3,-160/9]);
220 sub GSL_BLAS_DROTG : Tests {
221     my $a = [1];
222     my $b = [2];
224     my ($status, $c, $s) = gsl_blas_drotg($a, $b);
225     ok_similar( [$c, $s ], [ 1/sqrt(5) , 2/sqrt(5)  ] );
228 sub GSL_BLAS_DSYMV : Tests {
229     my $x = Math::GSL::Vector->new([1,2,3]);
230     my $y = Math::GSL::Vector->new([3,2,1]);
231     my $A = Math::GSL::Matrix->new(3,3);
232     map { gsl_matrix_set($A->raw, $_,0,$_+1); } (0..2);
233     gsl_matrix_set($A->raw, 0, 1, 2); 
234     gsl_matrix_set($A->raw, 1, 1, 1); 
235     gsl_matrix_set($A->raw, 2, 1, 2); 
236     gsl_matrix_set($A->raw, 0, 2, 3); 
237     gsl_matrix_set($A->raw, 1, 2, 2); 
238     gsl_matrix_set($A->raw, 2, 2, 1); 
239     ok_status(gsl_blas_dsymv($CblasLower, 2, $A->raw, $x->raw, 3, $y->raw));
240     ok_similar( [$y->as_list], [37,26,23]);
243 sub GSL_BLAS_DSYR : Tests {
244     my $x = Math::GSL::Vector->new([1,2,3]);
245     my $A = Math::GSL::Matrix->new(3,3);
246     gsl_matrix_set_zero($A->raw);
248     ok_status(gsl_blas_dsyr($CblasLower, 2, $x->raw, $A->raw));
249     ok_similar([ $A->row(0)->as_list ], [2,0,0]);
250     ok_similar([ $A->row(1)->as_list ], [4,8,0]);
251     ok_similar([ $A->row(2)->as_list ], [6,12,18]);
254 sub GSL_BLAS_ZHER : Tests {
255     my $x = gsl_vector_complex_alloc(2);
256     my $A = gsl_matrix_complex_alloc(2,2);
257     my $alpha = gsl_complex_rect(0,0);
259     for my $line (0 .. 1) {
260         map { gsl_matrix_complex_set($A, $_, $line, $alpha) } (0..1); 
261     }  
262     $alpha = gsl_complex_rect(1,2);
263     gsl_vector_complex_set($x, 0, $alpha);
264     $alpha = gsl_complex_rect(2,2);
265     gsl_vector_complex_set($x, 1, $alpha);
266     $alpha = gsl_complex_rect(1,1);
267     ok_status(gsl_blas_zher($CblasLower, 2, $x, $A));
269     ok_similar([ gsl_parts(gsl_matrix_complex_get($A,0,0)) ], [10, 0 ]);
270     ok_similar([ gsl_parts(gsl_matrix_complex_get($A,1,0)) ], [12, -4]);
271     ok_similar([ gsl_parts(gsl_matrix_complex_get($A,1,1)) ], [16, 0 ]);
272     ok_similar([ gsl_parts(gsl_matrix_complex_get($A,0,1)) ], [0 , 0 ]);
275 sub GSL_BLAS_DSYR2 : Tests {
276     my $x = Math::GSL::Vector->new([1,2,3]);
277     my $y = Math::GSL::Vector->new([3,2,1]);
278     my $A = Math::GSL::Matrix->new(3,3);
279     gsl_matrix_set_zero($A->raw);
280     map { gsl_matrix_set($A->raw, $_, 0, ($_+1)**2) } (0..2);
281     map { gsl_matrix_set($A->raw, $_, 1, ($_+1)**2) } (1..2);
282     gsl_matrix_set($A->raw, 0, 1, (1)**2);
283     gsl_matrix_set($A->raw, 0, 1, 4);
284     gsl_matrix_set($A->raw, 0, 2, 9);
285     gsl_matrix_set($A->raw, 1, 2, 9);
286     gsl_matrix_set($A->raw, 2, 2, 3);
287     ok_status(gsl_blas_dsyr2($CblasLower, 2, $x->raw, $y->raw, $A->raw));
288     my @got = $A->row(0)->as_list;
289     ok_similar([@got], [13, 4, 9]);
290     @got = $A->row(1)->as_list;
291     ok_similar([@got], [20, 20, 9]);
292     @got = $A->row(2)->as_list;
293     ok_similar([@got], [29, 25, 15]);
296 sub GSL_BLAS_DGEMM : Tests {
297     my $A = Math::GSL::Matrix->new(2,2);
298     gsl_matrix_set($A->raw, 0,0,1);
299     gsl_matrix_set($A->raw, 1,0,3);
300     gsl_matrix_set($A->raw, 0,1,4);
301     gsl_matrix_set($A->raw, 1,1,2);
302     my $B = Math::GSL::Matrix->new(2,2);
303     gsl_matrix_set($B->raw, 0,0,2);
304     gsl_matrix_set($B->raw, 1,0,5);
305     gsl_matrix_set($B->raw, 0,1,1);
306     gsl_matrix_set($B->raw, 1,1,3);
307     my $C = Math::GSL::Matrix->new(2,2);
308     gsl_matrix_set_zero($C->raw);
309     ok_status(gsl_blas_dgemm($CblasNoTrans, $CblasNoTrans, 1, $A->raw, $B->raw, 1, $C->raw));
310     my @got = $C->row(0)->as_list;
311     ok_similar([@got], [22, 13]);
312     @got = $C->row(1)->as_list;
313     ok_similar([@got], [16, 9]);
316 sub GSL_BLAS_DSYMM : Tests {
317     my $A = Math::GSL::Matrix->new(2,2);
318     gsl_matrix_set($A->raw, 0,0,1);
319     gsl_matrix_set($A->raw, 1,0,3);
320     gsl_matrix_set($A->raw, 0,1,4);
321     gsl_matrix_set($A->raw, 1,1,2);
322     my $B = Math::GSL::Matrix->new(2,2);
323     gsl_matrix_set($B->raw, 0,0,2);
324     gsl_matrix_set($B->raw, 1,0,5);
325     gsl_matrix_set($B->raw, 0,1,1);
326     gsl_matrix_set($B->raw, 1,1,3);
327     my $C = Math::GSL::Matrix->new(2,2);
328     gsl_matrix_set_zero($C->raw);
329     ok_status(gsl_blas_dsymm($CblasLeft, $CblasUpper, 1, $A->raw, $B->raw, 1, $C->raw));
330     my @got = $C->row(0)->as_list;
331     ok_similar([@got], [22, 13]);
332     @got = $C->row(1)->as_list;
333     ok_similar([@got], [18, 10]);
336 sub GSL_BLAS_ZGEMM : Tests {
337     my $A = gsl_matrix_complex_alloc(2,2);
338     my $alpha = gsl_complex_rect(1,2);
339     gsl_matrix_complex_set($A, 0,0,$alpha);
340     $alpha = gsl_complex_rect(3,1);
341     gsl_matrix_complex_set($A, 0,1,$alpha);
342     $alpha = gsl_complex_rect(2,2);
343     gsl_matrix_complex_set($A, 1,0,$alpha);
344     $alpha = gsl_complex_rect(0,2);
345     gsl_matrix_complex_set($A, 1,1,$alpha);
347     my $B = gsl_matrix_complex_alloc(2,2);
348     $alpha = gsl_complex_rect(1,1);
349     gsl_matrix_complex_set($B, 0,0,$alpha);
350     $alpha = gsl_complex_rect(2,2);
351     gsl_matrix_complex_set($B, 0,1,$alpha);
352     $alpha = gsl_complex_rect(1,2);
353     gsl_matrix_complex_set($B, 1,0,$alpha);
354     $alpha = gsl_complex_rect(1,3);
355     gsl_matrix_complex_set($B, 1,1,$alpha);
357     my $C = gsl_matrix_complex_alloc(2,2);
358     $alpha = gsl_complex_rect(0,0);
359     map { gsl_matrix_complex_set($C, 0, $_, $alpha) } (0..1);
360     map { gsl_matrix_complex_set($C, 1, $_, $alpha) } (0..1);
362     $alpha = gsl_complex_rect(2,0);
363     my $beta = gsl_complex_rect(1,0);
364     ok_status(gsl_blas_zgemm($CblasNoTrans, $CblasNoTrans, $alpha, $A, $B, $beta, $C));
365     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 0, 0))], [0,20]);
366     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 0, 1))], [-4,32]);
367     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 1, 0))], [-8,12]);
368     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 1, 1))], [-12,20]);
371 sub GSL_BLAS_ZSYMM : Tests {
372     my $A = gsl_matrix_complex_alloc(2,2);
373     my $alpha = gsl_complex_rect(1,2);
374     gsl_matrix_complex_set($A, 0,0,$alpha);
375     $alpha = gsl_complex_rect(3,1);
376     gsl_matrix_complex_set($A, 0,1,$alpha);
377     $alpha = gsl_complex_rect(3,1);
378     gsl_matrix_complex_set($A, 1,0,$alpha);
379     $alpha = gsl_complex_rect(0,2);
380     gsl_matrix_complex_set($A, 1,1,$alpha);
382     my $B = gsl_matrix_complex_alloc(2,2);
383     $alpha = gsl_complex_rect(1,1);
384     gsl_matrix_complex_set($B, 0,0,$alpha);
385     $alpha = gsl_complex_rect(2,2);
386     gsl_matrix_complex_set($B, 0,1,$alpha);
387     $alpha = gsl_complex_rect(2,2);
388     gsl_matrix_complex_set($B, 1,0,$alpha);
389     $alpha = gsl_complex_rect(1,3);
390     gsl_matrix_complex_set($B, 1,1,$alpha);
392     my $C = gsl_matrix_complex_alloc(2,2);
393     $alpha = gsl_complex_rect(0,0);
394     map { gsl_matrix_complex_set($C, 0, $_, $alpha) } (0..1);
395     map { gsl_matrix_complex_set($C, 1, $_, $alpha) } (0..1);
396     
397     $alpha = gsl_complex_rect(2,0);
398     my $beta = gsl_complex_rect(1,0);
399     ok_status(gsl_blas_zsymm($CblasLeft, $CblasUpper, $alpha, $A, $B, $beta, $C));
400     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 0, 0))], [6,22]);
401     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 0, 1))], [-4,32]);
402     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 1, 0))], [-4,16]);
403     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 1, 1))], [-4,20]);
405 sub GSL_BLAS_ZHEMM : Tests {
406     my $A = gsl_matrix_complex_alloc(2,2);
407     my $alpha = gsl_complex_rect(3,0);
408     gsl_matrix_complex_set($A, 0,0,$alpha);
409     $alpha = gsl_complex_rect(2,1);
410     gsl_matrix_complex_set($A, 0,1,$alpha);
411     $alpha = gsl_complex_rect(2,-1);
412     gsl_matrix_complex_set($A, 1,0,$alpha);
413     $alpha = gsl_complex_rect(1,0);
414     gsl_matrix_complex_set($A, 1,1,$alpha);
416     my $B = gsl_matrix_complex_alloc(2,2);
417     $alpha = gsl_complex_rect(1,0);
418     gsl_matrix_complex_set($B, 0,0,$alpha);
419     $alpha = gsl_complex_rect(2,2);
420     gsl_matrix_complex_set($B, 0,1,$alpha);
421     $alpha = gsl_complex_rect(2,-2);
422     gsl_matrix_complex_set($B, 1,0,$alpha);
423     $alpha = gsl_complex_rect(2,0);
424     gsl_matrix_complex_set($B, 1,1,$alpha);
426     my $C = gsl_matrix_complex_alloc(2,2);
427     $alpha = gsl_complex_rect(0,0);
428     map { gsl_matrix_complex_set($C, 0, $_, $alpha) } (0..1);
429     map { gsl_matrix_complex_set($C, 1, $_, $alpha) } (0..1);
430     
431     $alpha = gsl_complex_rect(2,0);
432     my $beta = gsl_complex_rect(1,0);
433     ok_status(gsl_blas_zhemm($CblasLeft, $CblasUpper, $alpha, $A, $B, $beta, $C));
434     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 0, 0))], [18,-4]);
435     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 0, 1))], [20,16]);
436     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 1, 0))], [8,-6]);
437     ok_similar([ gsl_parts(gsl_matrix_complex_get($C, 1, 1))], [16,4]);
440 sub GSL_BLAS_DTRMM : Tests {
441     my $A = Math::GSL::Matrix->new(2,2);
442     gsl_matrix_set($A->raw, 0,0,1);
443     gsl_matrix_set($A->raw, 1,0,3);
444     gsl_matrix_set($A->raw, 0,1,4);
445     gsl_matrix_set($A->raw, 1,1,2);
446     my $B = Math::GSL::Matrix->new(2,2);
447     gsl_matrix_set($B->raw, 0,0,2);
448     gsl_matrix_set($B->raw, 1,0,5);
449     gsl_matrix_set($B->raw, 0,1,1);
450     gsl_matrix_set($B->raw, 1,1,3);
451     my $C = Math::GSL::Matrix->new(2,2);
452     gsl_matrix_set_zero($C->raw);
453     ok_status(gsl_blas_dtrmm($CblasLeft, $CblasUpper, $CblasNoTrans, $CblasUnit, 1, $A->raw, $B->raw));
454     my @got = $B->row(0)->as_list;
455     ok_similar([@got], [22, 13]);
456     @got = $B->row(1)->as_list;
457     ok_similar([@got], [5, 3]);
460 sub GSL_BLAS_ZTRMM : Tests {
461     my $A = gsl_matrix_complex_alloc(2,2);
462     my $alpha = gsl_complex_rect(3,1);
463     gsl_matrix_complex_set($A, 0,1,$alpha);
465     my $B = gsl_matrix_complex_alloc(2,2);
466     $alpha = gsl_complex_rect(1,0);
467     gsl_matrix_complex_set($B, 0,0,$alpha);
468     $alpha = gsl_complex_rect(2,2);
469     gsl_matrix_complex_set($B, 0,1,$alpha);
470     $alpha = gsl_complex_rect(1,-2);
471     gsl_matrix_complex_set($B, 1,0,$alpha);
472     $alpha = gsl_complex_rect(4,2);
473     gsl_matrix_complex_set($B, 1,1,$alpha);
474     
475     $alpha = gsl_complex_rect(1,0);
476     ok_status(gsl_blas_ztrmm($CblasLeft, $CblasUpper, $CblasNoTrans, $CblasUnit, $alpha, $A, $B));
477     ok_similar([gsl_parts(gsl_matrix_complex_get($B, 0, 0))], [6, -5]);
478     ok_similar([gsl_parts(gsl_matrix_complex_get($B, 0, 1))], [12, 12]);
479     ok_similar([gsl_parts(gsl_matrix_complex_get($B, 1, 0))], [1, -2]);
480     ok_similar([gsl_parts(gsl_matrix_complex_get($B, 1, 1))], [4, 2]);
483 sub GSL_BLAS_DSYRK : Tests {
484     my $A = Math::GSL::Matrix->new(2,2);
485     gsl_matrix_set($A->raw, 0, 0, 1);
486     gsl_matrix_set($A->raw, 0, 1, 4);
487     gsl_matrix_set($A->raw, 1, 0, 4);
488     gsl_matrix_set($A->raw, 1, 1, 3);
489     my $C = Math::GSL::Matrix->new(2,2);
490     gsl_matrix_set_zero($C->raw);
491     ok_status(gsl_blas_dsyrk ($CblasUpper, $CblasNoTrans, 1, $A->raw, 1, $C->raw));
492     ok_similar([$C->row(0)->as_list], [17,16]);
493     ok_similar([$C->row(1)->as_list], [0,25]);
496 sub GSL_BLAS_ZSYRK : Tests {
497     my $A = gsl_matrix_complex_alloc(2,2);
498     my $alpha = gsl_complex_rect(3,1);
499     gsl_matrix_complex_set($A, 0,0,$alpha);
500     $alpha = gsl_complex_rect(2,1);
501     gsl_matrix_complex_set($A, 0,1,$alpha);
502     gsl_matrix_complex_set($A, 1,0,$alpha);
503     $alpha = gsl_complex_rect(1,1);
504     gsl_matrix_complex_set($A, 1,1,$alpha);
506     my $C = gsl_matrix_complex_alloc(2,2);
507     $alpha = gsl_complex_rect(0,0);
508     map { gsl_matrix_complex_set($C, 0,$_,$alpha) } (0..1);
509     map { gsl_matrix_complex_set($C, 1,$_,$alpha) } (0..1);
510     
511     $alpha = gsl_complex_rect(1,0);
512     my $beta = gsl_complex_rect(1,0);
513     ok_status(gsl_blas_zsyrk($CblasUpper, $CblasNoTrans, $alpha, $A, $beta, $C));
514     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 0, 0))], [11, 10]);
515     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 0, 1))], [6, 8]);
516     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 1, 0))], [0, 0]);
517     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 1, 1))], [3, 6]);
520 sub GSL_BLAS_ZHERK : Tests {
521     my $A = gsl_matrix_complex_alloc(2,2);
522     my $alpha = gsl_complex_rect(3,0);
523     gsl_matrix_complex_set($A, 0,0,$alpha);
524     $alpha = gsl_complex_rect(2,1);
525     gsl_matrix_complex_set($A, 0,1,$alpha);
526     $alpha = gsl_complex_rect(2,-1);
527     gsl_matrix_complex_set($A, 1,0,$alpha);
528     $alpha = gsl_complex_rect(1,0);
529     gsl_matrix_complex_set($A, 1,1,$alpha);
531     my $C = gsl_matrix_complex_alloc(2,2);
532     $alpha = gsl_complex_rect(0,0);
533     map { gsl_matrix_complex_set($C, 0,$_,$alpha) } (0..1);
534     map { gsl_matrix_complex_set($C, 1,$_,$alpha) } (0..1);
536     ok_status(gsl_blas_zherk ($CblasUpper, $CblasNoTrans, 1, $A, 1, $C));
537     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 0, 0))], [14, 0]);
538     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 0, 1))], [8, 4]);
539     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 1, 0))], [0, 0]);
540     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 1, 1))], [6, 0]);
543 sub GSL_BLAS_ZHER2K : Tests {
544     my $A = gsl_matrix_complex_alloc(2,2);
545     my $alpha = gsl_complex_rect(3,0);
546     gsl_matrix_complex_set($A, 0,0,$alpha);
547     $alpha = gsl_complex_rect(2,1);
548     gsl_matrix_complex_set($A, 0,1,$alpha);
549     $alpha = gsl_complex_rect(2,-1);
550     gsl_matrix_complex_set($A, 1,0,$alpha);
551     $alpha = gsl_complex_rect(1,0);
552     gsl_matrix_complex_set($A, 1,1,$alpha);
554     my $B = gsl_matrix_complex_alloc(2,2);
555     $alpha = gsl_complex_rect(6,0);
556     gsl_matrix_complex_set($A, 0,0,$alpha);
557     $alpha = gsl_complex_rect(3,1);
558     gsl_matrix_complex_set($A, 0,1,$alpha);
559     $alpha = gsl_complex_rect(3,-1);
560     gsl_matrix_complex_set($A, 1,0,$alpha);
561     $alpha = gsl_complex_rect(5,0);
562     gsl_matrix_complex_set($A, 1,1,$alpha);
563     
564     my $C = gsl_matrix_complex_alloc(2,2);
565     $alpha = gsl_complex_rect(0,0);
566     map { gsl_matrix_complex_set($C, 0,$_,$alpha) } (0..1);
567     map { gsl_matrix_complex_set($C, 1,$_,$alpha) } (0..1);
569     $alpha = gsl_complex_rect(1,0);
571     ok_status(gsl_blas_zher2k($CblasUpper, $CblasNoTrans, $alpha, $A, $B, 1, $C));
572     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 1, 0))], [0, 0]);
573     local $TODO = "These results follow the formula given by the documentation, don't know why it fails";
574     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 0, 0))], [50, 0]);
575     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 0, 1))], [34, 15]);
576     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 1, 1))], [24, 0]);
579 sub GSL_BLAS_DSYR2K : Tests {
580     my $A = Math::GSL::Matrix->new(2,2);
581     gsl_matrix_set($A->raw, 0, 0, 1);
582     gsl_matrix_set($A->raw, 0, 1, 4);
583     gsl_matrix_set($A->raw, 1, 0, 4);
584     gsl_matrix_set($A->raw, 1, 1, 3);
585     my $B = Math::GSL::Matrix->new(2,2);
586     gsl_matrix_set($B->raw, 0, 0, 2);
587     gsl_matrix_set($B->raw, 0, 1, 5);
588     gsl_matrix_set($B->raw, 1, 0, 5);
589     gsl_matrix_set($B->raw, 1, 1, 1);
590     my $C = Math::GSL::Matrix->new(2,2);
591     gsl_matrix_set_zero($C->raw);
592     ok_status(gsl_blas_dsyr2k ($CblasUpper, $CblasNoTrans, 1, $A->raw, $B->raw, 1, $C->raw, ));
593     ok_similar([$C->row(0)->as_list], [44,32]);
594     ok_similar([$C->row(1)->as_list], [0,46]);
597 sub GSL_BLAS_ZSYR2K : Tests {
598     my $A = gsl_matrix_complex_alloc(2,2);
599     my $alpha = gsl_complex_rect(3,0);
600     gsl_matrix_complex_set($A, 0,0,$alpha);
601     $alpha = gsl_complex_rect(2,1);
602     gsl_matrix_complex_set($A, 0,1,$alpha);
603     $alpha = gsl_complex_rect(2,1);
604     gsl_matrix_complex_set($A, 1,0,$alpha);
605     $alpha = gsl_complex_rect(1,0);
606     gsl_matrix_complex_set($A, 1,1,$alpha);
608     my $B = gsl_matrix_complex_alloc(2,2);
609     $alpha = gsl_complex_rect(6,0);
610     gsl_matrix_complex_set($A, 0,0,$alpha);
611     $alpha = gsl_complex_rect(3,1);
612     gsl_matrix_complex_set($A, 0,1,$alpha);
613     $alpha = gsl_complex_rect(3,1);
614     gsl_matrix_complex_set($A, 1,0,$alpha);
615     $alpha = gsl_complex_rect(5,0);
616     gsl_matrix_complex_set($A, 1,1,$alpha);
617     
618     my $C = gsl_matrix_complex_alloc(2,2);
619     $alpha = gsl_complex_rect(0,0);
620     map { gsl_matrix_complex_set($C, 0,$_,$alpha) } (0..1);
621     map { gsl_matrix_complex_set($C, 1,$_,$alpha) } (0..1);
623     $alpha = gsl_complex_rect(1,0);
624     my $beta = gsl_complex_rect(1,0);
626     ok_status(gsl_blas_zsyr2k($CblasUpper, $CblasNoTrans, $alpha, $A, $B, $beta, $C));
627     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 1, 0))], [0, 0]);
628     local $TODO = "These results follow the formula given by the documentation, don't know why it fails";
629     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 0, 0))], [46, 10]);
630     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 0, 1))], [34, 15]);
631     ok_similar([gsl_parts(gsl_matrix_complex_get($C, 1, 1))], [24, 4]);
635 Test::Class->runtests;