6 use Scalar::Util 'blessed';
8 use Math::GSL qw/:all/;
9 use Math::GSL::Errno qw/:all/;
10 use Math::GSL::Eigen qw/:all/;
11 use Math::GSL::Vector qw/:all/;
12 use Math::GSL::Complex qw/:all/;
13 use Math::GSL::BLAS qw/gsl_blas_dgemm/;
14 use Math::GSL::CBLAS qw/$CblasNoTrans/;
15 use Math::GSL::Permutation;
16 use Math::GSL::Linalg qw/
23 # should only include needed methods
24 use Math::GSL::MatrixComplex qw/:all/;
25 use Math::GSL::VectorComplex qw/:all/;
29 '*' => \&_multiplication,
35 gsl_matrix_alloc gsl_matrix_calloc gsl_matrix_alloc_from_block
36 gsl_matrix_alloc_from_matrix gsl_vector_alloc_row_from_matrix
37 gsl_vector_alloc_col_from_matrix gsl_matrix_free gsl_matrix_submatrix
38 gsl_matrix_row gsl_matrix_column gsl_matrix_diagonal
39 gsl_matrix_subdiagonal gsl_matrix_superdiagonal gsl_matrix_subrow
40 gsl_matrix_subcolumn gsl_matrix_view_array
41 gsl_matrix_view_array_with_tda gsl_matrix_view_vector
42 gsl_matrix_view_vector_with_tda gsl_matrix_const_submatrix
43 gsl_matrix_const_row gsl_matrix_const_column gsl_matrix_const_diagonal
44 gsl_matrix_const_subdiagonal gsl_matrix_const_superdiagonal
45 gsl_matrix_const_subrow gsl_matrix_const_subcolumn
46 gsl_matrix_const_view_array gsl_matrix_const_view_array_with_tda
47 gsl_matrix_const_view_vector gsl_matrix_const_view_vector_with_tda
48 gsl_matrix_get gsl_matrix_set gsl_matrix_ptr gsl_matrix_const_ptr
49 gsl_matrix_set_zero gsl_matrix_set_identity gsl_matrix_set_all
50 gsl_matrix_fread gsl_matrix_fwrite gsl_matrix_fscanf gsl_matrix_fprintf
51 gsl_matrix_memcpy gsl_matrix_swap gsl_matrix_swap_rows
52 gsl_matrix_swap_columns gsl_matrix_swap_rowcol gsl_matrix_transpose
53 gsl_matrix_transpose_memcpy gsl_matrix_max gsl_matrix_min
54 gsl_matrix_minmax gsl_matrix_max_index gsl_matrix_min_index
55 gsl_matrix_minmax_index gsl_matrix_isnull gsl_matrix_ispos
56 gsl_matrix_isneg gsl_matrix_isnonneg gsl_matrix_add gsl_matrix_sub
57 gsl_matrix_mul_elements gsl_matrix_div_elements gsl_matrix_scale
58 gsl_matrix_add_constant gsl_matrix_add_diagonal
59 gsl_matrix_char_alloc gsl_matrix_char_calloc gsl_matrix_char_alloc_from_block
60 gsl_matrix_char_alloc_from_matrix gsl_vector_char_alloc_row_from_matrix gsl_vector_char_alloc_col_from_matrix
61 gsl_matrix_char_free gsl_matrix_char_submatrix
62 gsl_matrix_char_row gsl_matrix_char_column
63 gsl_matrix_char_diagonal gsl_matrix_char_subdiagonal gsl_matrix_char_superdiagonal
64 gsl_matrix_char_subrow gsl_matrix_char_subcolumn gsl_matrix_char_view_array
65 gsl_matrix_char_view_array_with_tda gsl_matrix_char_view_vector gsl_matrix_char_view_vector_with_tda
66 gsl_matrix_char_const_submatrix gsl_matrix_char_const_row gsl_matrix_char_const_column
67 gsl_matrix_char_const_diagonal gsl_matrix_char_const_subdiagonal gsl_matrix_char_const_superdiagonal
68 gsl_matrix_char_const_subrow gsl_matrix_char_const_subcolumn gsl_matrix_char_const_view_array
69 gsl_matrix_char_const_view_array_with_tda gsl_matrix_char_const_view_vector gsl_matrix_char_const_view_vector_with_tda
70 gsl_matrix_char_get gsl_matrix_char_set gsl_matrix_char_ptr gsl_matrix_char_const_ptr
71 gsl_matrix_char_set_zero gsl_matrix_char_set_identity
72 gsl_matrix_char_set_all gsl_matrix_char_fread
73 gsl_matrix_char_fwrite gsl_matrix_char_fscanf gsl_matrix_char_fprintf
74 gsl_matrix_char_memcpy gsl_matrix_char_swap
75 gsl_matrix_char_swap_rows gsl_matrix_char_swap_columns
76 gsl_matrix_char_swap_rowcol gsl_matrix_char_transpose gsl_matrix_char_transpose_memcpy
77 gsl_matrix_char_max gsl_matrix_char_min
78 gsl_matrix_char_minmax gsl_matrix_char_max_index
79 gsl_matrix_char_min_index gsl_matrix_char_minmax_index
80 gsl_matrix_char_isnull gsl_matrix_char_ispos gsl_matrix_char_isneg
81 gsl_matrix_char_isnonneg gsl_matrix_char_add
82 gsl_matrix_char_sub gsl_matrix_char_mul_elements gsl_matrix_char_div_elements
83 gsl_matrix_char_scale gsl_matrix_char_add_constant gsl_matrix_char_add_diagonal
84 gsl_matrix_int_alloc gsl_matrix_int_calloc gsl_matrix_int_alloc_from_block
85 gsl_matrix_int_alloc_from_matrix gsl_vector_int_alloc_row_from_matrix gsl_vector_int_alloc_col_from_matrix
86 gsl_matrix_int_free gsl_matrix_int_submatrix gsl_matrix_int_row
87 gsl_matrix_int_column gsl_matrix_int_diagonal gsl_matrix_int_subdiagonal
88 gsl_matrix_int_superdiagonal gsl_matrix_int_subrow gsl_matrix_int_subcolumn gsl_matrix_int_view_array
89 gsl_matrix_int_view_array_with_tda gsl_matrix_int_view_vector gsl_matrix_int_view_vector_with_tda
90 gsl_matrix_int_const_submatrix gsl_matrix_int_const_row gsl_matrix_int_const_column
91 gsl_matrix_int_const_diagonal gsl_matrix_int_const_subdiagonal gsl_matrix_int_const_superdiagonal
92 gsl_matrix_int_const_subrow gsl_matrix_int_const_subcolumn gsl_matrix_int_const_view_array
93 gsl_matrix_int_const_view_array_with_tda gsl_matrix_int_const_view_vector gsl_matrix_int_const_view_vector_with_tda
94 gsl_matrix_int_get gsl_matrix_int_set
95 gsl_matrix_int_ptr gsl_matrix_int_const_ptr
96 gsl_matrix_int_set_zero gsl_matrix_int_set_identity gsl_matrix_int_set_all
97 gsl_matrix_int_fread gsl_matrix_int_fwrite
98 gsl_matrix_int_fscanf gsl_matrix_int_fprintf
99 gsl_matrix_int_memcpy gsl_matrix_int_swap
100 gsl_matrix_int_swap_rows gsl_matrix_int_swap_columns gsl_matrix_int_swap_rowcol
101 gsl_matrix_int_transpose gsl_matrix_int_transpose_memcpy
102 gsl_matrix_int_max gsl_matrix_int_min gsl_matrix_int_minmax
103 gsl_matrix_int_max_index gsl_matrix_int_min_index
104 gsl_matrix_int_minmax_index gsl_matrix_int_isnull
105 gsl_matrix_int_ispos gsl_matrix_int_isneg gsl_matrix_int_isnonneg
106 gsl_matrix_int_add gsl_matrix_int_sub
107 gsl_matrix_int_mul_elements gsl_matrix_int_div_elements gsl_matrix_int_scale
108 gsl_matrix_int_add_constant gsl_matrix_int_add_diagonal
109 gsl_matrix_get_row gsl_matrix_get_col gsl_matrix_set_row gsl_matrix_set_col
112 our %EXPORT_TAGS = ( all => [ @EXPORT_OK ],
114 gsl_matrix_char_alloc
115 gsl_matrix_char_calloc
116 gsl_matrix_char_alloc_from_block
117 gsl_matrix_char_alloc_from_matrix
118 gsl_vector_char_alloc_row_from_matrix
119 gsl_vector_char_alloc_col_from_matrix
121 gsl_matrix_char_submatrix
123 gsl_matrix_char_column
124 gsl_matrix_char_diagonal
125 gsl_matrix_char_subdiagonal
126 gsl_matrix_char_superdiagonal
127 gsl_matrix_char_subrow
128 gsl_matrix_char_subcolumn
129 gsl_matrix_char_view_array
130 gsl_matrix_char_view_array_with_tda
131 gsl_matrix_char_view_vector
132 gsl_matrix_char_view_vector_with_tda
133 gsl_matrix_char_const_submatrix
134 gsl_matrix_char_const_row
135 gsl_matrix_char_const_column
136 gsl_matrix_char_const_diagonal
137 gsl_matrix_char_const_subdiagonal
138 gsl_matrix_char_const_superdiagonal
139 gsl_matrix_char_const_subrow
140 gsl_matrix_char_const_subcolumn
141 gsl_matrix_char_const_view_array
142 gsl_matrix_char_const_view_array_with_tda
143 gsl_matrix_char_const_view_vector
144 gsl_matrix_char_const_view_vector_with_tda
148 gsl_matrix_char_const_ptr
149 gsl_matrix_char_set_zero
150 gsl_matrix_char_set_identity
151 gsl_matrix_char_set_all
152 gsl_matrix_char_fread
153 gsl_matrix_char_fwrite
154 gsl_matrix_char_fscanf
155 gsl_matrix_char_fprintf
156 gsl_matrix_char_memcpy
158 gsl_matrix_char_swap_rows
159 gsl_matrix_char_swap_columns
160 gsl_matrix_char_swap_rowcol
161 gsl_matrix_char_transpose
162 gsl_matrix_char_transpose_memcpy
165 gsl_matrix_char_minmax
166 gsl_matrix_char_max_index
167 gsl_matrix_char_min_index
168 gsl_matrix_char_minmax_index
169 gsl_matrix_char_isnull
170 gsl_matrix_char_ispos
171 gsl_matrix_char_isneg
172 gsl_matrix_char_isnonneg
175 gsl_matrix_char_mul_elements
176 gsl_matrix_char_div_elements
177 gsl_matrix_char_scale
178 gsl_matrix_char_add_constant
179 gsl_matrix_char_add_diagonal
185 gsl_matrix_alloc_from_block
186 gsl_matrix_alloc_from_matrix
187 gsl_vector_alloc_row_from_matrix
188 gsl_vector_alloc_col_from_matrix
194 gsl_matrix_subdiagonal
195 gsl_matrix_superdiagonal
198 gsl_matrix_view_array
199 gsl_matrix_view_array_with_tda
200 gsl_matrix_view_vector
201 gsl_matrix_view_vector_with_tda
202 gsl_matrix_const_submatrix
204 gsl_matrix_const_column
205 gsl_matrix_const_diagonal
206 gsl_matrix_const_subdiagonal
207 gsl_matrix_const_superdiagonal
208 gsl_matrix_const_subrow
209 gsl_matrix_const_subcolumn
210 gsl_matrix_const_view_array
211 gsl_matrix_const_view_array_with_tda
212 gsl_matrix_const_view_vector
213 gsl_matrix_const_view_vector_with_tda
219 gsl_matrix_set_identity
228 gsl_matrix_swap_columns
229 gsl_matrix_swap_rowcol
231 gsl_matrix_transpose_memcpy
236 gsl_matrix_minmax_index
242 gsl_matrix_mul_elements
243 gsl_matrix_div_elements
245 gsl_matrix_add_constant
246 gsl_matrix_add_diagonal
250 gsl_matrix_int_alloc_from_matrix
252 gsl_matrix_int_column
253 gsl_matrix_int_superdiagonal
254 gsl_matrix_int_view_array_with_tda
255 gsl_matrix_int_const_submatrix
256 gsl_matrix_int_const_diagonal
257 gsl_matrix_int_const_subrow
258 gsl_matrix_int_const_view_array_with_tda
261 gsl_matrix_int_set_zero
263 gsl_matrix_int_fscanf
264 gsl_matrix_int_memcpy
265 gsl_matrix_int_swap_rows
266 gsl_matrix_int_transpose
268 gsl_matrix_int_max_index
269 gsl_matrix_int_minmax_index
272 gsl_matrix_int_mul_elements
273 gsl_matrix_int_add_constant
280 Math::GSL::Matrix - Mathematical functions concerning Matrices
284 use Math::GSL::Matrix qw/:all/;
285 my $matrix1 = Math::GSL::Matrix->new(5,5); # OO interface
286 my $matrix2 = $matrix1 + 4; # You can add or substract values or matrices to OO matrices
287 my $matrix3 = $matrix1 - 4;
288 my $matrix4 = $matrix2 + $matrix1;
289 my $matrix5 = $matrix2 . $matrix1; # This is a scalar product, it simply multiply each element
290 # with the element of $matrix1 that have the same position
291 # See Math::GSL::BLAS if you want scalar product
293 my $matrix6 = $matrix2 . 8; # Multiply every elements of $matrix2 by 8
294 my $matrix7 = $matrix2 * $matrix1; # scalar product of two matrices
295 my $matrix8 = gsl_matrix_alloc(5,5); # standard interface
298 =head1 Objected Oriented Interface to GSL Math::GSL::Matrix
300 =head2 Math::GSL::Matrix->new()
302 Creates a new Matrix of the given size.
304 my $matrix = Math::GSL::Matrix->new(10,10);
309 my ($class, $rows, $cols) = @_;
312 if ( defined $rows && defined $cols &&
313 $rows > 0 && $cols > 0 &&
314 (int $rows == $rows) && (int $cols == $cols)){
316 $matrix = gsl_matrix_alloc($rows,$cols);
318 croak( __PACKAGE__.'::new($x,$y) - $x and $y must be positive integers');
320 gsl_matrix_set_zero($matrix);
321 $this->{_matrix} = $matrix;
322 ($this->{_rows}, $this->{_cols}) = ($rows,$cols);
327 Get the underlying GSL matrix object created by SWIG, useful for using gsl_matrix_* functions which do not have an OO counterpart.
329 my $matrix = Math::GSL::Matrix->new(3,3);
330 my $gsl_matrix = $matrix->raw;
331 my $stuff = gsl_matrix_get($gsl_matrix, 1, 2);
334 sub raw { (shift)->{_matrix} }
338 Set a matrix to the identity matrix, i.e. one on the diagonal and zero elsewhere.
340 my $I = $matrix->identity;
347 gsl_matrix_set_identity($self->raw);
353 Set a matrix to the zero matrix.
362 gsl_matrix_set_zero($self->raw);
368 Returns a copy of the matrix, which has the same size and values but resides at a different location in memory.
370 my $matrix = Math::GSL::Matrix->new(5,5);
371 my $copy = $matrix->copy;
377 my $copy = Math::GSL::Matrix->new( $self->rows, $self->cols );
378 if ( gsl_matrix_memcpy($copy->raw, $self->raw) != $GSL_SUCCESS ) {
379 croak "Math::GSL - error copying memory, aborting";
387 Returns the number of rows in the matrix.
389 my $rows = $matrix->rows;
392 sub rows { (shift)->{_rows} }
396 Returns the number of columns in the matrix.
398 my $cols = $matrix->cols;
401 sub cols { (shift)->{_cols} }
405 Get the contents of a Math::GSL::Matrix object as a Perl list.
407 my $matrix = Math::GSL::Matrix->new(3,3);
409 my @matrix = $matrix->as_list;
415 Returns true if a matrix is square, i.e. it has the same number of rows as columns, false otherwise.
422 return ($self->rows == $self->cols) ? 1 : 0 ;
427 Returns the determinant of a matrix (computed by LU decomposition) or dies if called on a non-square matrix.
429 my $det = $matrix->det();
436 croak(__PACKAGE__."- determinant only exists for square matrices") unless $self->is_square;
438 my $p = Math::GSL::Permutation->new( $self->rows );
439 my $LU = $self->copy;
441 my ($result, $s) = gsl_linalg_LU_decomp($LU->raw, $p->raw);
442 return gsl_linalg_LU_det($LU->raw, $s );
447 Returns the natural log of the absolute value of the determinant of a matrix (computed by LU decomposition) or dies if called on a non-square matrix.
449 my $lndet = $matrix->lndet();
456 croak(__PACKAGE__."- log determinant only exists for square matrices") unless $self->is_square;
458 my $p = Math::GSL::Permutation->new( $self->rows );
459 my $LU = $self->copy;
461 gsl_linalg_LU_decomp($LU->raw, $p->raw);
462 return gsl_linalg_LU_lndet($LU->raw);
467 Returns the inverse of a matrix or dies when called on a non-square matrix.
469 my $inverse = $matrix->inverse;
476 croak(__PACKAGE__."- inverse only exists for square matrices") unless $self->is_square;
478 my $p = Math::GSL::Permutation->new( $self->rows );
479 my $LU = $self->copy;
480 my $inverse = $self->copy;
482 # should check return status
483 gsl_linalg_LU_decomp($LU->raw, $p->raw);
484 gsl_linalg_LU_invert($LU->raw, $p->raw,$inverse->raw);
498 my ($r,$c) = ($self->rows,$self->cols);
499 croak "Math::GSL::Matrix : \$matrix->eigenvalues - \$matrix must be square" unless ($r == $c);
500 my $evec = gsl_matrix_complex_alloc($r,$c);
501 my $eigen = gsl_eigen_nonsymmv_alloc($r);
502 my $vector = gsl_vector_complex_alloc($r);
504 gsl_eigen_nonsymmv($self->raw, $vector, $evec, $eigen);
506 my $x = gsl_vector_complex_real($vector); # vector of real components
507 my $y = gsl_vector_complex_imag($vector); # vector of imaginary components
510 _as_complex( gsl_vector_get($x->{vector}, $_) ,
511 gsl_vector_get($y->{vector}, $_) )
519 my ($r,$c) = ($self->rows,$self->cols);
520 croak "Math::GSL::Matrix : \$matrix->eigenvalues - \$matrix must be square" unless ($r == $c);
521 my $evec = Math::GSL::MatrixComplex->new($r,$c);
522 my $eigen = gsl_eigen_nonsymmv_alloc($r);
523 my $vector = Math::GSL::VectorComplex->new($r);
525 gsl_eigen_nonsymmv($self->raw, $vector->raw, $evec->raw, $eigen);
527 my $eigenvectors = [ map { $evec->col($_)->as_vector } (0,1) ];
529 my $x = gsl_vector_complex_real($vector->raw); # vector of real components
530 my $y = gsl_vector_complex_imag($vector->raw); # vector of imaginary components
533 _as_complex( gsl_vector_get($x->{vector}, $_) ,
534 gsl_vector_get($y->{vector}, $_) )
537 return ($eigenvalues, $eigenvectors);
542 my ($x,$y) = ($w,$v);
543 if( ref $w eq 'Math::GSL::Complex') {
544 ($x,$y) = (gsl_real($w),gsl_imag($w));
546 my $z = Math::Complex->make( $x, $y);
551 Returns a 1xN or Nx1 matrix as a Math::GSL::Vector object. Dies if called on a matrix that is not a single row or column. Useful for turning the output of C<col()> or C<row()> into a vector, like so:
553 my $vector1 = $matrix->col(0)->as_vector;
554 my $vector2 = $matrix->row(1)->as_vector;
561 croak(__PACKAGE__.'::as_vector() - must be a single row or column matrix') unless ($self->cols == 1 or $self->rows == 1);
563 # TODO: there is a faster way to do this
564 return Math::GSL::Vector->new([ $self->as_list ] );
570 Returns a matrix as a flattened Perl list.
572 my @values = $matrix->as_list;
581 for my $row ( 0 .. $self->rows-1) {
582 push @matrix, map { gsl_matrix_get($self->raw, $row, $_) } (0 .. $self->cols-1 );
589 Returns a row matrix of the row you enter.
591 my $matrix = Math::GSL::Matrix->new(3,3);
593 my $matrix_row = $matrix->row(0);
599 my ($self, $row) = @_;
600 croak (__PACKAGE__.'::$matrix->row($row) - invalid $row value')
601 unless (($row < $self->rows) and $row >= 0);
603 my $rowmat = Math::GSL::Matrix->new(1,$self->cols);
605 for my $n (0 .. $self->cols-1) {
606 gsl_matrix_set( $rowmat->raw, 0, $n,
607 gsl_matrix_get($self->raw, $row, $n)
616 Returns a col matrix of the column you enter.
618 my $matrix = Math::GSL::Matrix->new(3,3);
620 my $matrix_col = $matrix->col(0);
626 my ($self, $col) = @_;
627 croak (__PACKAGE__."::\$matrix->col(\$col) - $col not a valid column")
628 unless ($col < $self->cols and $col >= 0);
630 my $colmat = Math::GSL::Matrix->new($self->rows, 1);
632 map { gsl_matrix_set($colmat->raw, $_, 0,
633 gsl_matrix_get($self->raw, $_, $col),
635 } (0 .. $self->rows-1);
642 Sets a the values of a row with the elements of an array.
644 my $matrix = Math::GSL::Matrix->new(3,3);
645 $matrix->set_row(0, [8, 6, 2]);
647 You can also set multiple rows at once with chained calls:
648 my $matrix = Math::GSL::Matrix->new(3,3);
649 $matrix->set_row(0, [8, 6, 2])
650 ->set_row(1, [2, 4, 1]);
656 my ($self, $row, $values) = @_;
657 my $length = $#$values;
658 die __PACKAGE__.'::new($x, $values) - $values must be a nonempty array reference' if $length == -1;
659 die __PACKAGE__.'::set_row($x, $values) - $x must be a valid row number'
660 if ($row < 0 || $row >= $self->rows);
661 die __PACKAGE__.'::set_row($x, $values) - $values must contains the same number of elements as there is columns in the matrix'
662 if($length != $self->cols-1);
663 map { gsl_matrix_set($self->raw, $row, $_, $values->[$_]) } (0..$length);
669 Sets a the values of a column with the elements of an array.
671 my $matrix = Math::GSL::Matrix->new(3,3);
672 $matrix->set_col(0, [8, 6, 2]);
674 You can also set multiple columns at once with chained calls:
675 my $matrix = Math::GSL::Matrix->new(3,3);
676 $matrix->set_col(0, [8, 6, 2])
677 ->set_col(1, [2, 4, 1]);
683 my ($self, $col, $values) = @_;
684 my $length = $#$values;
686 die __PACKAGE__.'::new($x, $values) - $values must be a nonempty array reference' if $length == -1;
687 die __PACKAGE__.'::set_row($x, $values) - $x must be a valid column number'
688 if ($col < 0 || $col >= $self->cols);
689 die __PACKAGE__.'::set_row($x, $values) - $values must contains the same number of elements as there is rowss in the matrix'
690 if($length != $self->rows-1);
692 map { gsl_matrix_set($self->raw, $_, $col, $values->[$_]) } (0..$length);
697 my ($left, $right) = @_;
699 my $lcopy = $left->copy;
701 if ( blessed $right && $right->isa('Math::GSL::Matrix') && blessed $left && $left->isa('Math::GSL::Matrix') ) {
702 if ( $left->cols == $right->cols && $left->rows == $right->rows ) {
703 gsl_matrix_add($lcopy->raw, $right->raw);
705 croak "Math::GSL - addition of matrices must be called with two objects matrices and must have the same number of columns and rows";
708 gsl_matrix_add_constant($lcopy->raw, $right);
714 my ($left, $right) = @_;
716 my $lcopy = $left->copy;
718 if ( blessed $right && $right->isa('Math::GSL::Matrix') && blessed $left && $left->isa('Math::GSL::Matrix') ) {
719 if ( $left->cols == $right->cols && $left->rows == $right->rows ) {
720 gsl_matrix_sub($lcopy->raw, $right->raw);
722 croak "Math::GSL - subtraction of matrices must be called with two objects matrices and must have the same number of columns and rows";
725 gsl_matrix_add_constant($lcopy->raw, $right*-1);
730 sub _multiplication {
731 my ($left,$right) = @_;
732 my $lcopy = $left->copy;
734 if ( blessed $right && $right->isa('Math::GSL::Matrix') ) {
735 if ( $left->cols == $right->cols && $left->rows == $right->rows ) {
736 return _dot_product($left,$right);
737 #gsl_matrix_mul_elements($lcopy->raw, $right->raw);
739 croak "Math::GSL::Matrix - multiplication of elements of matrices must be called with two objects matrices and must have the same number of columns and rows";
742 gsl_matrix_scale($lcopy->raw, $right);
748 my ($left,$right) = @_;
750 croak "Math::GSL::Matrix - incompatible matrices in multiplication"
751 unless ( blessed $right && $right->isa('Math::GSL::Matrix') and $left->rows == $right->cols );
752 my $C = Math::GSL::Matrix->new($left->rows, $right->cols);
753 gsl_blas_dgemm($CblasNoTrans, $CblasNoTrans, 1, $left->raw, $right->raw, 1, $C->raw);
759 Here is a list of all the functions included in this module :
763 =item C<gsl_matrix_alloc($i, $j)> - Return a gsl_matrix of $i rows and $j columns
765 =item C<gsl_matrix_calloc($i, $j)> - Return a gsl_matrix of $i rows and $j columns and initialize all of the elements of the matrix to zero
767 =item C<gsl_matrix_alloc_from_block> -
769 =item C<gsl_matrix_free> -
771 =item C<gsl_matrix_alloc_from_matrix > -
773 =item C<gsl_vector_alloc_row_from_matrix> -
775 =item C<gsl_vector_alloc_col_from_matrix > -
777 =item C<gsl_matrix_submatrix($m, $k1, $k2, $n1, $n2)> - Return a matrix view of the matrix $m. The upper-left element of the submatrix is the element ($k1,$k2) of the original matrix. The submatrix has $n1 rows and $n2 columns.
779 =item C<gsl_matrix_row($m , $i)> - Return a vector view of the $i-th row of the matrix $m
781 =item C<gsl_matrix_column($m, $j)> - Return a vector view of the $j-th column of the matrix $m
783 =item C<gsl_matrix_diagonal($m)> - Return a vector view of the diagonal of the vector. The matrix doesn't have to be square.
785 =item C<gsl_matrix_subdiagonal($m, $k)> - Return a vector view of the $k-th subdiagonal of the matrix $m. The diagonal of the matrix corresponds to k=0.
787 =item C<gsl_matrix_superdiagonal($m, $k)> - Return a vector view of the $k-th superdiagonal of the matrix $m. The matrix doesn't have to be square.
789 =item C<gsl_matrix_subrow($m, $i, $offset, $n)> - Return a vector view of the $i-th row of the matrix $m beginning at offset elements and containing n elements.
791 =item C<gsl_matrix_subcolumn($m, $j, $offset, $n)> - Return a vector view of the $j-th column of the matrix $m beginning at offset elements and containing n elements.
793 =item C<gsl_matrix_view_array($base, $n1, $n2)> - This function returns a matrix view of the array reference $base. The matrix has $n1 rows and $n2 columns. The physical number of columns in memory is also given by $n2. Mathematically, the (i,j)-th element of the new matrix is given by, m'(i,j) = $base->[i*$n2 + j] where the index i runs from 0 to $n1-1 and the index j runs from 0 to $n2-1. The new matrix is only a view of the array reference $base. When the view goes out of scope the original array reference $base will continue to exist. The original memory can only be deallocated by freeing the original array. Of course, the original array should not be deallocated while the view is still in use.
795 =item C<gsl_matrix_view_array_with_tda($base, $n1, $n2, $tda)> - This function returns a matrix view of the array reference $base with a physical number of columns $tda which may differ from the corresponding dimension of the matrix. The matrix has $n1 rows and $n2 columns, and the physical number of columns in memory is given by $tda. Mathematically, the (i,j)-th element of the new matrix is given by, m'(i,j) = $base->[i*$tda + j] where the index i runs from 0 to $n1-1 and the index j runs from 0 to $n2-1. The new matrix is only a view of the array reference $base. When the view goes out of scope the original array reference $base will continue to exist. The original memory can only be deallocated by freeing the original array. Of course, the original array should not be deallocated while the view is still in use.
797 =item C<gsl_matrix_view_vector> -
799 =item C<gsl_matrix_view_vector_with_tda> -
801 =item C<gsl_matrix_const_submatrix> -
803 =item C<gsl_matrix_get($m, $i, $j)> - Return the (i,j)-th element of the matrix $m
805 =item C<gsl_matrix_set($m, $i, $j, $x)> - Set the value of the (i,j)-th element of the matrix $m to $x
807 =item C<gsl_matrix_ptr> -
809 =item C<gsl_matrix_const_ptr> -
811 =item C<gsl_matrix_set_zero($m)> - Set all the elements of the matrix $m to zero
813 =item C<gsl_matrix_set_identity($m)> - Set the elements of the matrix $m to the corresponding elements of the identity matrix
815 =item C<gsl_matrix_set_all($m, $x)> - Set all the elements of the matrix $m to the value $x
817 =item C<gsl_matrix_fread($fh, $m)> - Read a file which has been written with gsl_matrix_fwrite from the stream $fh opened with the gsl_fopen function from the Math::GSL module and stores the data inside the matrix $m
819 =item C<gsl_matrix_fwrite($fh, $m)> - Write the elements of the matrix $m in binary format to a stream $fh opened with the gsl_fopen function from the Math::GSL module
821 =item C<gsl_matrix_fscanf($fh, $m)> - Read a file which has been written with gsl_matrix_fprintf from the stream $fh opened with the gsl_fopenfunction from the Math::GSL module and stores the data inside the matrix $m
823 =item C<gsl_matrix_fprintf($fh, $m, $format)> - Write the elements of the matrix $m in the format $format (for example "%f" is the format for double) to a stream $fh opened with the gsl_fopen function from the Math::GSL module
825 =item C<gsl_matrix_memcpy($dest, $src)> - Copy the elements of the matrix $src to the matrix $dest. The two matrices must have the same size.
827 =item C<gsl_matrix_swap($m1, $m2)> - Exchange the elements of the matrices $m1 and $m2 by copying. The two matrices must have the same size.
829 =item C<gsl_matrix_swap_rows($m, $i, $j)> - Exchange the $i-th and $j-th row of the matrix $m. The function returns 0 if the operation suceeded, 1 otherwise.
831 =item C<gsl_matrix_swap_columns($m, $i, $j)> - Exchange the $i-th and $j-th column of the matrix $m. The function returns 0 if the operation suceeded, 1 otherwise.
833 =item C<gsl_matrix_swap_rowcol($m, $i, $j)> - Exchange the $i-th row and the $j-th column of the matrix $m. The matrix must be square. The function returns 0 if the operation suceeded, 1 otherwise.
835 =item C<gsl_matrix_transpose($m)> - This function replaces the matrix m by its transpose by copying the elements of the matrix in-place. The matrix must be square for this operation to be possible.
837 =item C<gsl_matrix_transpose_memcpy($dest, $src)> - Make the matrix $dest the transpose of the matrix $src. This function works for all matrices provided that the dimensions of the matrix dest match the transposed dimensions of the matrix src.
839 =item C<gsl_matrix_max($m)> - Return the maximum value in the matrix $m
841 =item C<gsl_matrix_min($m)> - Return the minimum value in the matrix $m
843 =item C<gsl_matrix_minmax($m)> - Return two values, the first is the minimum value of the Matrix $m and the second is the maximum of the same the same matrix.
845 =item C<gsl_matrix_max_index($m)> - Return two values, the first is the the i indice of the maximum value of the matrix $m and the second is the j indice of the same value.
847 =item C<gsl_matrix_min_index($m)> - Return two values, the first is the the i indice of the minimum value of the matrix $m and the second is the j indice of the same value.
849 =item C<gsl_matrix_minmax_index($m)> - Return four values, the first is the i indice of the minimum of the matrix $m, the second is the j indice of the same value, the third is the i indice of the maximum of the matrix $m and the fourth is the j indice of the same value
851 =item C<gsl_matrix_isnull($m)> - Return 1 if all the elements of the matrix $m are zero, 0 otherwise
853 =item C<gsl_matrix_ispos($m)> - Return 1 if all the elements of the matrix $m are strictly positve, 0 otherwise
855 =item C<gsl_matrix_isneg($m)> - Return 1 if all the elements of the matrix $m are strictly negative, 0 otherwise
857 =item C<gsl_matrix_isnonneg($m)> - Return 1 if all the elements of the matrix $m are non-negatuive, 0 otherwise
859 =item C<gsl_matrix_add($a, $b)> - Add the elements of matrix $b to the elements of matrix $a
861 =item C<gsl_matrix_sub($a, $b)> - Subtract the elements of matrix $b from the elements of matrix $a
863 =item C<gsl_matrix_mul_elements($a, $b)> - Multiplie the elements of matrix $a by the elements of matrix $b
865 =item C<gsl_matrix_div_elements($a, $b)> - Divide the elements of matrix $a by the elements of matrix $b
867 =item C<gsl_matrix_scale($a, $x)> - Multiplie the elements of matrix $a by the constant factor $x
869 =item C<gsl_matrix_add_constant($a, $x)> - Add the constant value $x to the elements of the matrix $a
871 =item C<gsl_matrix_add_diagonal($a, $x)> - Add the constant value $x to the elements of the diagonal of the matrix $a
873 =item C<gsl_matrix_get_row($v, $m, $i)> - Copy the elements of the $i-th row of the matrix $m into the vector $v. The lenght of the vector must be of the same as the lenght of the row. The function returns 0 if it succeded, 1 otherwise.
875 =item C<gsl_matrix_get_col($v, $m, $i)> - Copy the elements of the $j-th column of the matrix $m into the vector $v. The lenght of the vector must be of the same as the lenght of the column. The function returns 0 if it succeded, 1 otherwise.
877 =item C<gsl_matrix_set_row($m, $i, $v)> - Copy the elements of vector $v into the $i-th row of the matrix $m The lenght of the vector must be of the same as the lenght of the row. The function returns 0 if it succeded, 1 otherwise.
879 =item C<gsl_matrix_set_col($m, $j, $v)> - Copy the elements of vector $v into the $j-th row of the matrix $m The lenght of the vector must be of the same as the lenght of the column. The function returns 0 if it succeded, 1 otherwise.
883 These are related to constant views of a matrix.
888 =item C<gsl_matrix_const_row>
890 =item C<gsl_matrix_const_colum>
892 =item C<gsl_matrix_const_diagonal>
894 =item C<gsl_matrix_const_subdiagonal>
896 =item C<gsl_matrix_const_superdiagonal>
898 =item C<gsl_matrix_const_subrow>
900 =item C<gsl_matrix_const_subcolumn>
902 =item C<gsl_matrix_const_view_array>
904 =item C<gsl_matrix_const_view_array_with_tda>
909 The following functions are similar to those above but work with C<char>'s and C<int>'s. We are not quite
910 sure if anyone wants these. Please speak up if you do and/or submit some patches to this documentation, please!
915 =item gsl_matrix_const_view_vector
917 =item gsl_matrix_const_view_vector_with_tda
919 =item gsl_matrix_char_alloc
921 =item gsl_matrix_char_calloc
923 =item gsl_matrix_char_alloc_from_block
925 =item gsl_matrix_char_alloc_from_matrix
927 =item gsl_vector_char_alloc_row_from_matrix
929 =item gsl_vector_char_alloc_col_from_matrix
931 =item gsl_matrix_char_free
933 =item gsl_matrix_char_submatrix
935 =item gsl_matrix_char_row
937 =item gsl_matrix_char_column
939 =item gsl_matrix_char_diagonal
941 =item gsl_matrix_char_subdiagonal
943 =item gsl_matrix_char_superdiagonal
945 =item gsl_matrix_char_subrow
947 =item gsl_matrix_char_subcolumn
949 =item gsl_matrix_char_view_array
951 =item gsl_matrix_char_view_array_with_tda
953 =item gsl_matrix_char_view_vector
955 =item gsl_matrix_char_view_vector_with_tda
957 =item gsl_matrix_char_const_submatrix
959 =item gsl_matrix_char_const_row
961 =item gsl_matrix_char_const_column
963 =item gsl_matrix_char_const_diagonal
965 =item gsl_matrix_char_const_subdiagonal
967 =item gsl_matrix_char_const_superdiagonal
969 =item gsl_matrix_char_const_subrow
971 =item gsl_matrix_char_const_subcolumn
973 =item gsl_matrix_char_const_view_array
975 =item gsl_matrix_char_const_view_array_with_tda
977 =item gsl_matrix_char_const_view_vector
979 =item gsl_matrix_char_const_view_vector_with_tda
981 =item gsl_matrix_char_get
983 =item gsl_matrix_char_set
985 =item gsl_matrix_char_ptr
987 =item gsl_matrix_char_const_ptr
989 =item gsl_matrix_char_set_zero
991 =item gsl_matrix_char_set_identity
993 =item gsl_matrix_char_set_all
995 =item gsl_matrix_char_fread
997 =item gsl_matrix_char_fwrite
999 =item gsl_matrix_char_fscanf
1001 =item gsl_matrix_char_fprintf
1003 =item gsl_matrix_char_memcpy
1005 =item gsl_matrix_char_swap
1007 =item gsl_matrix_char_swap_rows
1009 =item gsl_matrix_char_swap_columns
1011 =item gsl_matrix_char_swap_rowcol
1013 =item gsl_matrix_char_transpose
1015 =item gsl_matrix_char_transpose_memcpy
1017 =item gsl_matrix_char_max
1019 =item gsl_matrix_char_min
1021 =item gsl_matrix_char_minmax
1023 =item gsl_matrix_char_max_index
1025 =item gsl_matrix_char_min_index
1027 =item gsl_matrix_char_minmax_index
1029 =item gsl_matrix_char_isnull
1031 =item gsl_matrix_char_ispos
1033 =item gsl_matrix_char_isneg
1035 =item gsl_matrix_char_isnonneg
1037 =item gsl_matrix_char_add
1039 =item gsl_matrix_char_sub
1041 =item gsl_matrix_char_mul_elements
1043 =item gsl_matrix_char_div_elements
1045 =item gsl_matrix_char_scale
1047 =item gsl_matrix_char_add_constant
1049 =item gsl_matrix_char_add_diagonal
1051 =item gsl_matrix_int_alloc
1053 =item gsl_matrix_int_calloc
1055 =item gsl_matrix_int_alloc_from_block
1057 =item gsl_matrix_int_alloc_from_matrix
1059 =item gsl_vector_int_alloc_row_from_matrix
1061 =item gsl_vector_int_alloc_col_from_matrix
1063 =item gsl_matrix_int_free
1065 =item gsl_matrix_int_submatrix
1067 =item gsl_matrix_int_row
1069 =item gsl_matrix_int_column
1071 =item gsl_matrix_int_diagonal
1073 =item gsl_matrix_int_subdiagonal
1075 =item gsl_matrix_int_superdiagonal
1077 =item gsl_matrix_int_subrow
1079 =item gsl_matrix_int_subcolumn
1081 =item gsl_matrix_int_view_array
1083 =item gsl_matrix_int_view_array_with_tda
1085 =item gsl_matrix_int_view_vector
1087 =item gsl_matrix_int_view_vector_with_tda
1089 =item gsl_matrix_int_const_submatrix
1091 =item gsl_matrix_int_const_row
1093 =item gsl_matrix_int_const_column
1095 =item gsl_matrix_int_ptr
1097 =item gsl_matrix_int_const_ptr
1099 =item gsl_matrix_int_set_zero
1101 =item gsl_matrix_int_set_identity
1103 =item gsl_matrix_int_set_all
1105 =item gsl_matrix_int_fread
1107 =item gsl_matrix_int_fwrite
1109 =item gsl_matrix_int_fscanf
1111 =item gsl_matrix_int_fprintf
1113 =item gsl_matrix_int_memcpy
1115 =item gsl_matrix_int_swap
1117 =item gsl_matrix_int_swap_rows
1119 =item gsl_matrix_int_swap_columns
1121 =item gsl_matrix_int_swap_rowcol
1123 =item gsl_matrix_int_transpose
1125 =item gsl_matrix_int_transpose_memcpy
1127 =item gsl_matrix_int_max
1129 =item gsl_matrix_int_min
1131 =item gsl_matrix_int_minmax
1133 =item gsl_matrix_int_max_index
1135 =item gsl_matrix_int_min_index
1137 =item gsl_matrix_int_minmax_index
1139 =item gsl_matrix_int_isnull
1141 =item gsl_matrix_int_ispos
1143 =item gsl_matrix_int_isneg
1145 =item gsl_matrix_int_isnonneg
1147 =item gsl_matrix_int_add
1149 =item gsl_matrix_int_sub
1151 =item gsl_matrix_int_mul_elements
1153 =item gsl_matrix_int_div_elements
1155 =item gsl_matrix_int_scale
1157 =item gsl_matrix_int_add_constant
1159 =item gsl_matrix_int_add_diagonal
1163 You have to add the functions you want to use inside the qw /put_funtion_here /.
1164 You can also write use Math::GSL::Matrix qw/:all/ to use all avaible functions of the module.
1165 Other tags are also avaible, here is a complete list of all tags for this module :
1181 For more informations on the functions, we refer you to the GSL offcial documentation
1182 L<http://www.gnu.org/software/gsl/manual/html_node/>
1184 Tip : search on google: site:http://www.gnu.org/software/gsl/manual/html_node/name_of_the_function_you_want
1188 Most of the examples from this section are perl versions of the examples at L<http://www.gnu.org/software/gsl/manual/html_node/Example-programs-for-matrices.html>
1190 The program below shows how to allocate, initialize and read from a matrix using the functions gsl_matrix_alloc, gsl_matrix_set and gsl_matrix_get.
1192 use Math::GSL::Matrix qw/:all/;
1193 my $m = gsl_matrix_alloc (10,3);
1196 gsl_matrix_set($m, $i, $j, 0.23 + 100*$i + $j);
1200 for my $i (0..99){ # OUT OF RANGE ERROR
1202 print "m($i, $j) = " . gsl_matrix_get ($m, $i, $j) . "\n";
1205 gsl_matrix_free ($m);
1208 use Math::GSL::Matrix qw/:all/;
1210 my $m = gsl_matrix_alloc (100, 100);
1211 my $a = gsl_matrix_alloc (100, 100);
1215 gsl_matrix_set ($m, $i, $j, 0.23 + $i + $j);
1219 The next program shows how to write a matrix to a file.
1221 my $out = gsl_fopen("test.dat", "wb");
1222 gsl_matrix_fwrite ($out, $m);
1225 my $in = gsl_fopen("test.dat", "rb");
1226 gsl_matrix_fread ($in, $a);
1232 $mij = gsl_matrix_get ($m, $i, $j);
1233 $aij = gsl_matrix_get ($a, $i, $j);
1234 $k++ if ($mij != $aij);
1238 gsl_matrix_free($m);
1239 gsl_matrix_free($a);
1241 print "differences = $k (should be zero)\n";
1245 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
1247 =head1 COPYRIGHT AND LICENSE
1249 Copyright (C) 2008 Jonathan Leto and Thierry Moisan
1251 This program is free software; you can redistribute it and/or modify it
1252 under the same terms as Perl itself.