Update META.yml
[Math-GSL.git] / lib / Math / GSL / Matrix / Test.pm
blob3c5a824ad1e1de1c8aefdeffdcce3d2e84ff67d4
1 package Math::GSL::Matrix::Test;
2 use Math::GSL::Test qw/:all/;
3 use base q{Test::Class};
4 use Test::More;
5 use Math::GSL::Matrix qw/:all/;
6 use Math::GSL::Vector qw/:all/;
7 use Math::GSL::Complex qw/:all/;
8 use Math::GSL::Errno qw/:all/;
9 use Math::GSL qw/:all/;
10 use Data::Dumper;
11 use strict;
13 BEGIN{ gsl_set_error_handler_off(); }
15 sub make_fixture : Test(setup) {
16 my $self = shift;
17 $self->{matrix} = gsl_matrix_alloc(5,5);
18 $self->{obj} = Math::GSL::Matrix->new(5,5);
19 gsl_matrix_set_zero($self->{matrix});
22 sub teardown : Test(teardown) {
23 unlink 'matrix' if -f 'matrix';
26 sub GSL_MATRIX_ALLOC : Tests {
27 my $matrix = gsl_matrix_alloc(5,5);
28 isa_ok($matrix, 'Math::GSL::Matrix');
31 sub GSL_MATRIX_SET : Tests {
32 my $self = shift;
33 map { gsl_matrix_set($self->{matrix}, $_,$_, $_ ** 2) } (0..4);
35 my @got = map { gsl_matrix_get($self->{matrix}, $_, $_) } (0..4);
36 ok_similar( [ @got ], [ map { $_**2 } (0..4) ] );
39 sub GSL_MATRIX_CALLOC : Tests {
40 my $matrix = gsl_matrix_calloc(5,5);
41 isa_ok($matrix, 'Math::GSL::Matrix');
43 my @got = map { gsl_matrix_get($matrix, $_, $_) } (0..4);
44 ok_similar( [ @got ], [ (0) x 5 ], 'gsl_matrix_calloc' );
47 sub GSL_MATRIX_FREE : Tests {
48 my $matrix = gsl_matrix_calloc(5,5);
49 isa_ok($matrix, 'Math::GSL::Matrix');
51 is(gsl_matrix_get($matrix, 0, 0), 0);
52 gsl_matrix_free($matrix);
55 sub GSL_MATRIX_SUBMATRIX : Tests {
56 my $matrix = gsl_matrix_alloc(5,5);
57 map { gsl_matrix_set($matrix, $_,$_, $_) } (0..4);
58 my $subMatrix = gsl_matrix_submatrix($matrix, 0, 0, 2, 2);
59 my @got = map { gsl_matrix_get($matrix, $_, $_) } (0..2);
60 ok_similar( [ @got ], [ 0..2 ] );
63 sub GSL_MATRIX_ROW : Tests {
64 my $self = shift;
65 for my $line (0..4) {
66 map { gsl_matrix_set($self->{matrix}, $_,$line, $_) } (0..4);
69 my $vector_view = gsl_matrix_row($self->{matrix}, 2);
70 my @got = map { gsl_vector_get($vector_view->{vector}, $_) } (0..4);
71 ok_similar( [ @got ], [ (2)x 5], 'gsl_matrix_row' );
74 sub GSL_MATRIX_COLUMN : Tests {
75 my $self = shift;
76 my $view = gsl_vector_alloc(5);
78 for my $line (0..4) {
79 map { gsl_matrix_set($self->{matrix}, $line,$_, $line*$_) } (0..4);
81 $view = gsl_matrix_column($self->{matrix}, 2);
82 my $vec = $view->swig_vector_get();
84 my @got = map { gsl_vector_get($vec, $_) } (0..4);
85 ok_similar( [ @got ], [0,2,4,6,8 ], 'gsl_matrix_column' );
88 sub GSL_MATRIX_DIAGONAL : Tests {
89 my $matrix = gsl_matrix_alloc(4,4);
90 map { gsl_matrix_set($matrix, $_,$_, $_) } (0..3);
91 my $view = gsl_matrix_diagonal($matrix);
92 # better interface is needed
93 my $vec = $view->swig_vector_get();
95 my @got = map { gsl_vector_get($vec, $_) } (0..3);
96 ok_similar( [ @got ], [ 0 .. 3 ], 'gsl_matrix_diagonal');
99 sub GSL_MATRIX_SUBDIAGONAL : Tests {
100 my $matrix = gsl_matrix_alloc(4,4);
102 map { gsl_matrix_set($matrix, $_,$_, $_) } (0..3);
104 my $view = gsl_matrix_subdiagonal($matrix, 0);
105 my $vec = $view->swig_vector_get();
106 my @got = map { gsl_vector_get($vec, $_) } (0..3);
107 ok_similar( [ @got ], [ 0 .. 3 ], 'gsl_matrix_subdiagonal');
110 sub GSL_MATRIX_SWAP : Tests {
111 my $self=shift;
112 map { gsl_matrix_set($self->{matrix}, $_,$_, $_ ** 2) } (0..4);
113 my $matrix = gsl_matrix_alloc(5,5);
114 map { gsl_matrix_set($matrix, $_,$_, $_) } (0..4);
115 ok_status(gsl_matrix_swap($self->{matrix}, $matrix));
117 my @got = map { gsl_matrix_get($self->{matrix}, $_, $_) } (0..4);
118 map { is($got[$_], $_) } (0..4);
119 @got = map { gsl_matrix_get($matrix, $_, $_) } (0..4);
120 map { is($got[$_], $_** 2) } (0..4);
123 sub GSL_MATRIX_MEMCPY : Tests {
124 my $self = shift;
125 my $matrix = gsl_matrix_alloc(5,5);
126 map { gsl_matrix_set($self->{matrix}, $_,$_, $_ ** 2) } (0..4);
127 ok_status(gsl_matrix_memcpy($matrix, $self->{matrix}));
128 ok_similar( [ map { gsl_matrix_get($matrix, $_, $_) } (0..4) ],
129 [ map { $_** 2 } (0..4) ]
133 sub GSL_MATRIX_SWAP_ROWS : Tests {
134 my $self = shift;
135 map { gsl_matrix_set($self->{matrix}, 0,$_, $_) } (0..4);
136 map { gsl_matrix_set($self->{matrix}, 1,$_, 3) } (0..4);
137 ok_status(gsl_matrix_swap_rows($self->{matrix}, 0, 1));
138 my @got = map { gsl_matrix_get($self->{matrix}, 1, $_) } (0..4);
139 map { is($got[$_], $_) } (0..4);
140 @got = map { gsl_matrix_get($self->{matrix}, 0, $_) } (0..4);
141 map { is($got[$_], 3) } (0..4);
144 sub GSL_MATRIX_SWAP_COLUMNS : Tests {
145 my $self = shift;
146 map { gsl_matrix_set($self->{matrix}, $_,0, $_) } (0..4);
147 map { gsl_matrix_set($self->{matrix}, $_,1, 3) } (0..4);
148 ok_status(gsl_matrix_swap_columns($self->{matrix}, 0, 1));
149 my @got = map { gsl_matrix_get($self->{matrix}, $_, 1) } (0..4);
150 map { is($got[$_], $_) } (0..4);
151 @got = map { gsl_matrix_get($self->{matrix}, $_, 0) } (0..4);
152 map { is($got[$_], 3) } (0..4);
155 sub GSL_MATRIX_SWAP_ROWCOL : Tests {
156 my $self = shift;
157 map { gsl_matrix_set($self->{matrix}, 0,$_, $_) } (0..4);
158 map { gsl_matrix_set($self->{matrix}, $_,2, 2) } (0..4);
159 ok_status(gsl_matrix_swap_rowcol($self->{matrix}, 0, 2));
161 my @got = map { gsl_matrix_get($self->{matrix}, $_, 2) } (0..4);
162 is_deeply( [ @got ], [ qw/2 1 0 3 4/ ] );
164 @got = map { gsl_matrix_get($self->{matrix}, 0, $_) } (0..4);
165 is_deeply( [ @got ], [ qw/2 2 2 2 2/ ] );
168 sub GSL_MATRIX_TRANSPOSE_MEMCPY : Tests {
169 my $self = shift;
170 map { gsl_matrix_set($self->{matrix}, 0,$_, $_) } (0..4);
171 my $matrix = gsl_matrix_alloc(5,5);
172 ok_status(gsl_matrix_transpose_memcpy($matrix, $self->{matrix}));
173 my @got = map { gsl_matrix_get($matrix, $_, 0) } (0..4);
174 map { is($got[$_], $_) } (0..4);
177 sub GSL_MATRIX_TRANSPOSE : Tests {
178 my $self = shift;
179 map { gsl_matrix_set($self->{matrix}, 0,$_, $_) } (0..4);
180 is(gsl_matrix_transpose($self->{matrix}), 0);
181 my @got = map { gsl_matrix_get($self->{matrix}, $_, 0) } (0..4);
182 map { is($got[$_], $_) } (0..4);
185 sub GSL_MATRIX_ADD : Tests {
186 my $self = shift;
187 my $matrix = gsl_matrix_alloc(5, 5);
188 map { gsl_matrix_set($self->{matrix}, $_, $_, $_) } (0..4);
189 map { gsl_matrix_set($matrix, $_, $_, $_) } (0..4);
190 is(gsl_matrix_add($self->{matrix}, $matrix), 0);
191 my @got = map { gsl_matrix_get($self->{matrix}, $_, $_) } (0..4);
192 map { is($got[$_], $_*2) } (0..4);
195 sub GSL_MATRIX_SUB : Tests {
196 my $self = shift;
197 my $matrix = gsl_matrix_alloc(5, 5);
198 map { gsl_matrix_set($self->{matrix}, $_, $_, $_) } (0..4);
199 map { gsl_matrix_set($matrix, $_, $_, $_) } (0..4);
200 is(gsl_matrix_sub($self->{matrix}, $matrix), 0);
201 my @got = map { gsl_matrix_get($self->{matrix}, $_, $_) } (0..4);
202 map { is($got[$_], 0) } (0..4);
205 sub GSL_MATRIX_MUL_ELEMENTS : Tests {
206 my $self = shift;
207 my $matrix = gsl_matrix_alloc(5, 5);
208 map { gsl_matrix_set($self->{matrix}, $_, $_, $_) } (0..4);
209 map { gsl_matrix_set($matrix, $_, $_, $_) } (0..4);
210 is(gsl_matrix_mul_elements($self->{matrix}, $matrix), 0);
211 my @got = map { gsl_matrix_get($self->{matrix}, $_, $_) } (0..4);
212 map { is($got[$_], $_**2) } (0..4);
215 sub GSL_MATRIX_DIV_ELEMENTS : Tests {
216 my $self = shift;
217 my $matrix = gsl_matrix_alloc(5, 5);
218 map { gsl_matrix_set($self->{matrix}, $_, $_, $_+1) } (0..4);
219 map { gsl_matrix_set($matrix, $_, $_, $_+1) } (0..4);
220 is(gsl_matrix_div_elements($self->{matrix}, $matrix), 0);
221 my @got = map { gsl_matrix_get($self->{matrix}, $_, $_) } (0..4);
222 map { is($got[$_], 1) } (0..4);
225 sub GSL_MATRIX_SCALE : Tests {
226 my $self = shift;
227 map { gsl_matrix_set($self->{matrix}, $_, $_, $_) } (0..4);
228 is(gsl_matrix_scale($self->{matrix}, 4), 0);
229 my @got = map { gsl_matrix_get($self->{matrix}, $_, $_) } (0..4);
230 map { is($got[$_], $_*4) } (0..4);
233 sub GSL_MATRIX_ADD_CONSTANT : Tests {
234 my $self = shift;
235 map { gsl_matrix_set($self->{matrix}, $_, $_, $_) } (0..4);
236 is(gsl_matrix_add_constant($self->{matrix}, 8), 0);
237 my @got = map { gsl_matrix_get($self->{matrix}, $_, $_) } (0..4);
238 map { is($got[$_], $_+8) } (0..4);
241 sub GSL_MATRIX_MAX : Tests {
242 my $self = shift;
243 for my $row (0..4) {
244 map { gsl_matrix_set($self->{matrix}, $row, $_, $_**2 ) } (0..4);
246 is(gsl_matrix_max($self->{matrix}), 16);
249 sub GSL_MATRIX_MIN : Tests {
250 my $self = shift;
251 for my $row (0..4) {
252 map { gsl_matrix_set($self->{matrix}, $row, $_, $_**2 ) } (0..4);
254 is(gsl_matrix_min($self->{matrix}), 0);
257 sub GSL_MATRIX_MINMAX : Test {
258 my $self = shift;
259 map { gsl_matrix_set($self->{matrix}, $_, $_, $_**2) } (0..4);
260 my ($min, $max) = gsl_matrix_minmax($self->{matrix});
261 ok_similar( [ $min, $max ], [ 0, 16], 'gsl_matrix_minmax' );
264 sub GSL_MATRIX_MAX_INDEX : Tests {
265 my $self = shift;
266 for my $row (0..3) {
267 map { gsl_matrix_set($self->{matrix}, $row, $_, $_) } (0..4);
269 map { gsl_matrix_set($self->{matrix}, $_, $_, $_**2) } (0..4);
270 my ($imax, $jmax) = gsl_matrix_max_index($self->{matrix});
271 ok_similar( [ $imax, $jmax ], [ 4, 4 ], 'gsl_matrix_max_index' );
274 sub GSL_MATRIX_MIN_INDEX : Tests {
275 my $self = shift;
276 map { gsl_matrix_set($self->{matrix}, $_, $_, $_**2) } (0..4);
277 my ($imin, $jmin) = gsl_matrix_min_index($self->{matrix});
278 ok_similar( [ $imin, $jmin ], [ 0, 0 ], 'gsl_matrix_min_index' );
282 sub GSL_MATRIX_ISNULL : Tests {
283 my $self = shift;
284 for my $line (0..4) {
285 map { gsl_matrix_set($self->{matrix}, $line, $_, 0) } (0..4);
287 is(gsl_matrix_isnull($self->{matrix}), 1);
289 for my $line (0..4) {
290 map { gsl_matrix_set($self->{matrix}, $line, $_, $_) } (0..4);
293 is(gsl_matrix_isnull($self->{matrix}), 0);
296 sub GSL_MATRIX_ISPOS : Tests {
297 my $self = shift;
298 my $line;
300 for($line=0; $line<5; $line++) {
301 map { gsl_matrix_set($self->{matrix}, $line, $_, -1) } (0..4); }
302 is(gsl_matrix_ispos($self->{matrix}), 0);
304 for($line=0; $line<5; $line++) {
305 map { gsl_matrix_set($self->{matrix}, $line, $_, 1) } (0..4); }
306 is(gsl_matrix_ispos($self->{matrix}), 1);
308 for($line=0; $line<5; $line++) {
309 map { gsl_matrix_set($self->{matrix}, $line, $_, 0) } (0..4); }
310 is(gsl_matrix_ispos($self->{matrix}), 0);
313 sub GSL_MATRIX_ISNEG : Tests {
314 my $self = shift;
315 my $line;
317 for($line=0; $line<5; $line++) {
318 map { gsl_matrix_set($self->{matrix}, $line, $_, -1) } (0..4); }
319 is(gsl_matrix_isneg($self->{matrix}), 1);
321 for($line=0; $line<5; $line++) {
322 map { gsl_matrix_set($self->{matrix}, $line, $_, 1) } (0..4); }
323 is(gsl_matrix_isneg($self->{matrix}), 0);
325 for($line=0; $line<5; $line++) {
326 map { gsl_matrix_set($self->{matrix}, $line, $_, 0) } (0..4); }
327 is(gsl_matrix_isneg($self->{matrix}), 0);
330 sub GSL_MATRIX_ISNONNEG : Tests {
331 my $self = shift;
332 for my $row (0..4) {
333 map { gsl_matrix_set($self->{matrix}, $row, $_, -1) } (0..4);
335 is(gsl_matrix_isnonneg($self->{matrix}), 0);
337 for my $row (0..4) {
338 map { gsl_matrix_set($self->{matrix}, $row, $_, 1) } (0..4);
340 is(gsl_matrix_isnonneg($self->{matrix}), 1);
342 for my $row (0..4) {
343 map { gsl_matrix_set($self->{matrix}, $row, $_, 0) } (0..4);
345 is(gsl_matrix_isnonneg($self->{matrix}), 1);
349 sub GSL_MATRIX_GET_ROW : Tests {
350 my $self = shift;
351 my $vector = gsl_vector_alloc(5);
352 map { gsl_matrix_set($self->{matrix}, 0, $_, $_) } (0..4);
353 is(gsl_matrix_get_row($vector, $self->{matrix}, 0), 0);
354 map { is(gsl_vector_get($vector, $_), $_) } (0..4);
357 sub GSL_MATRIX_GET_COL : Tests {
358 my $self = shift;
359 my $vector = gsl_vector_alloc(5);
360 map { gsl_matrix_set($self->{matrix}, $_, 0, $_) } (0..4);
361 is(gsl_matrix_get_col($vector, $self->{matrix}, 0), 0);
362 map { is(gsl_vector_get($vector, $_), $_) } (0..4);
365 sub GSL_MATRIX_SET_ROW : Tests {
366 my $self = shift;
367 my $vector = gsl_vector_alloc(5);
368 map { gsl_vector_set($vector, $_, $_**2) } (0..4);
369 is(gsl_matrix_set_row($self->{matrix}, 1, $vector), 0);
371 ok_similar( [ map { gsl_matrix_get($self->{matrix}, 1, $_) } (0..4) ],
372 [ map { $_ ** 2 } (0..4) ],
376 sub GSL_MATRIX_SET_COL : Tests {
377 my $self = shift;
378 my $vector = gsl_vector_alloc(5);
379 map { gsl_vector_set($vector, $_, $_**2) } (0..4);
380 is(gsl_matrix_set_col($self->{matrix}, 1, $vector), 0);
382 ok_similar( [ map { gsl_matrix_get($self->{matrix}, $_, 1) } (0..4) ],
383 [ map { $_ ** 2 } (0..4) ],
387 sub GSL_MATRIX_FREAD_FWRITE : Tests {
388 my $self = shift;
389 for my $line (0..4) {
390 map { gsl_matrix_set($self->{matrix}, $line, $_, $_**2) } (0..4);
393 my $fh = gsl_fopen("matrix", 'w');
394 is( gsl_matrix_fwrite($fh, $self->{matrix}), 0);
395 gsl_fclose($fh);
397 for my $line (0..4) {
398 map { gsl_matrix_set($self->{matrix}, $line, $_, $_**3) } (0..4);
401 $fh = gsl_fopen("matrix", 'r');
403 is(gsl_matrix_fread($fh, $self->{matrix}), 0);
404 for my $line (0..4) {
405 ok_similar(
406 [ map { gsl_matrix_get($self->{matrix}, $line, $_) } (0..4) ],
407 [ map { $_**2 } (0..4) ],
410 gsl_fclose($fh);
413 sub GSL_MATRIX_FPRINTF_FSCANF : Tests {
414 my $self = shift;
416 for my $line (0..4) {
417 map { gsl_matrix_set($self->{matrix}, $line, $_, $_**2) } (0..4);
420 my $fh = gsl_fopen("matrix", 'w');
421 is( gsl_matrix_fprintf($fh, $self->{matrix}, "%f"), 0);
422 ok_status(gsl_fclose($fh));
424 for my $line (0..4) {
425 map { gsl_matrix_set($self->{matrix}, $line, $_, $_**3) } (0..4);
428 $fh = gsl_fopen("matrix", 'r');
430 is(gsl_matrix_fscanf($fh, $self->{matrix}), 0);
431 for my $line (0..4) {
432 ok_similar(
433 [ map { gsl_matrix_get($self->{matrix}, $line, $_) } (0..4) ],
434 [ map { $_**2 } (0..4) ],
437 ok_status(gsl_fclose($fh));
440 sub GSL_MATRIX_MINMAX_INDEX : Tests {
441 my $self = shift;
442 my $line;
443 for ($line = 0; $line<4; $line ++) {
444 map { gsl_matrix_set($self->{matrix}, $line, $_, $_) } (0..4); }
445 map { gsl_matrix_set($self->{matrix}, 4, $_, $_**2) } (0..4);
446 my ($imin, $jmin, $imax, $jmax) = gsl_matrix_minmax_index($self->{matrix});
447 ok_similar( [ $imin, $jmin, $imax, $jmax ], [ 0, 0, 4, 4], 'gsl_matrix_minmax_index' );
450 sub GSL_MATRIX_ADD_DIAGONAL : Tests {
451 my $self = shift;
452 map { gsl_matrix_set($self->{matrix}, $_, $_, $_) } (0..4);
453 gsl_matrix_add_diagonal($self->{matrix}, 4);
454 ok_similar( [ map { gsl_matrix_get($self->{matrix}, $_, $_)} (0..4) ],
455 [ 4 .. 8 ],
459 sub GSL_MATRIX_NEW : Tests {
460 my $self = shift;
461 isa_ok( $self->{obj}, 'Math::GSL::Matrix' );
462 isa_ok( $self->{obj}->raw, 'Math::GSL::Matrix::gsl_matrix' );
463 ok( $self->{obj}->rows == 5, '->rows' );
464 ok( $self->{obj}->cols == 5, '->cols' );
467 sub AS_LIST_SQUARE : Tests {
468 my $matrix = Math::GSL::Matrix->new(5,5);
469 map { gsl_matrix_set($matrix->raw, $_, $_, 5 + $_**2) } (0..4);
470 is_deeply( [
471 5, 0, 0, 0, 0,
472 0, 6, 0, 0, 0,
473 0, 0, 9, 0, 0,
474 0, 0, 0,14, 0,
475 0, 0, 0, 0, 21
477 [ $matrix->as_list],
478 '$matrix->as_list',
482 sub AS_LIST_ROW : Tests {
483 my $matrix = Math::GSL::Matrix->new(1,5);
484 map { gsl_matrix_set($matrix->raw, 0 , $_, 5 + $_**2) } (0..4);
485 is_deeply( [ 5, 6, 9, 14, 21, ],
486 [ $matrix->as_list],
487 '$matrix->as_list',
491 sub ROW : Tests {
492 my $matrix = Math::GSL::Matrix->new(5,5);
493 map { gsl_matrix_set($matrix->raw, $_, $_, 5 + $_**2) } (0..4);
494 ok_similar( [qw/0 0 9 0 0/] , [$matrix->row(2)->as_list ] );
495 map { gsl_matrix_set($matrix->raw, 0, $_, $_) } (0..4);
496 ok_similar( [qw/0 1 2 3 4/] , [$matrix->row(0)->as_list ] );
500 sub COL : Tests {
501 my $matrix = Math::GSL::Matrix->new(3,3);
502 map { gsl_matrix_set($matrix->raw, $_, $_, 7 + $_**2) } (0..2);
503 ok_similar([7,0,0], [$matrix->col(0)->as_list]);
506 sub NEW_SETS_VALUES_TO_ZERO : Tests {
507 my $matrix = Math::GSL::Matrix->new(5,5);
508 my $sum;
510 map { $sum += $_ } $matrix->as_list;
511 ok( $sum == 0, 'new sets values to zero');
513 sub HERMITIAN : Tests {
514 my $matrix = gsl_matrix_complex_alloc(2,2);
515 my $transpose = gsl_matrix_complex_alloc(2,2);
516 gsl_matrix_complex_set($matrix, 0, 0, gsl_complex_rect(3,0));
517 gsl_matrix_complex_set($matrix, 0, 1, gsl_complex_rect(2,1));
518 gsl_matrix_complex_set($matrix, 1, 0, gsl_complex_rect(2,-1));
519 gsl_matrix_complex_set($matrix, 1, 1, gsl_complex_rect(1,0));
520 gsl_matrix_complex_memcpy($transpose, $matrix);
521 gsl_matrix_complex_transpose($transpose);
523 for my $row (0..1) {
524 map { gsl_matrix_complex_set($transpose, $row, $_, gsl_complex_conjugate(gsl_matrix_complex_get($transpose, $row, $_))) } (0..1);
527 my $upper_right = gsl_matrix_complex_get($matrix, 0, 1 );
528 my $lower_left = gsl_matrix_complex_get($matrix, 1, 0 );
530 ok( gsl_complex_eq( gsl_complex_conjugate($upper_right), $lower_left ), 'hermitian' );
533 sub SET_ROW : Tests {
534 my $m = Math::GSL::Matrix->new(3,3);
535 $m->set_row(0, [1,2,3]);
536 ok_similar([$m->row(0)->as_list], [1,2,3]);
538 sub SET_ROW_CHAINED : Tests {
539 my $m = Math::GSL::Matrix->new(3,3);
540 $m->set_row(1, [4,5,6])
541 ->set_row(2, [9,8,7]);
542 ok_similar([$m->row(1)->as_list], [4,5,6]);
543 ok_similar([$m->row(2)->as_list], [9,8,7]);
546 sub SET_COL_CHAINED : Tests {
547 my $m = Math::GSL::Matrix->new(3,3);
548 $m->set_col(1, [4,5,6])
549 ->set_col(2, [9,8,7]);
550 ok_similar([$m->col(1)->as_list], [4,5,6]);
551 ok_similar([$m->col(2)->as_list], [9,8,7]);
554 sub GSL_MATRIX_VIEW_ARRAY : Tests {
555 my $array = [8,4,3,7];
556 my $matrix_view = gsl_matrix_view_array ($array, 2,2);
557 ok_similar([map { gsl_matrix_get($matrix_view->{matrix}, 0, $_) } 0..1], [8, 4]);
558 ok_similar([map { gsl_matrix_get($matrix_view->{matrix}, 1, $_) } 0..1], [3, 7]);
561 sub GSL_MATRIX_VIEW_ARRAY_WITH_TDA : Tests {
562 my $array = [8,4,3,7,5];
563 my $matrix_view = gsl_matrix_view_array_with_tda ($array, 2,2, 3);
564 ok_similar([map { gsl_matrix_get($matrix_view->{matrix}, 0, $_) } 0..1], [8, 4]);
565 ok_similar([map { gsl_matrix_get($matrix_view->{matrix}, 1, $_) } 0..1], [7, 5]);