4 use Scalar::Util 'blessed';
5 use Math::GSL qw/:all/;
6 use Math::GSL::Complex qw/:all/;
7 use Math::GSL::Errno qw/:all/;
8 use Math::GSL::Permutation;
9 use Math::GSL::Eigen qw/:all/;
10 use Math::GSL::VectorComplex qw/:all/;
13 use Math::GSL::BLAS qw/gsl_blas_zgemm/;
14 use Math::GSL::CBLAS qw/$CblasNoTrans/;
15 use Math::GSL::Linalg qw/
16 gsl_linalg_complex_LU_decomp
17 gsl_linalg_complex_LU_det
18 gsl_linalg_complex_LU_lndet
19 gsl_linalg_complex_LU_invert
22 '*' => \&_multiplication,
32 gsl_matrix_complex_alloc
33 gsl_matrix_complex_calloc
34 gsl_matrix_complex_alloc_from_block
35 gsl_matrix_complex_alloc_from_matrix
36 gsl_vector_complex_alloc_row_from_matrix
37 gsl_vector_complex_alloc_col_from_matrix
38 gsl_matrix_complex_free
39 gsl_matrix_complex_submatrix
40 gsl_matrix_complex_row
41 gsl_matrix_complex_column
42 gsl_matrix_complex_diagonal
43 gsl_matrix_complex_subdiagonal
44 gsl_matrix_complex_superdiagonal
45 gsl_matrix_complex_subrow
46 gsl_matrix_complex_subcolumn
47 gsl_matrix_complex_view_array
48 gsl_matrix_complex_view_array_with_tda
49 gsl_matrix_complex_view_vector
50 gsl_matrix_complex_view_vector_with_tda
51 gsl_matrix_complex_const_submatrix
52 gsl_matrix_complex_const_row
53 gsl_matrix_complex_const_column
54 gsl_matrix_complex_const_diagonal
55 gsl_matrix_complex_const_subdiagonal
56 gsl_matrix_complex_const_superdiagonal
57 gsl_matrix_complex_const_subrow
58 gsl_matrix_complex_const_subcolumn
59 gsl_matrix_complex_const_view_array
60 gsl_matrix_complex_const_view_array_with_tda
61 gsl_matrix_complex_const_view_vector
62 gsl_matrix_complex_const_view_vector_with_tda
63 gsl_matrix_complex_get
64 gsl_matrix_complex_set
65 gsl_matrix_complex_ptr
66 gsl_matrix_complex_const_ptr
67 gsl_matrix_complex_set_zero
68 gsl_matrix_complex_set_identity
69 gsl_matrix_complex_set_all
70 gsl_matrix_complex_fread
71 gsl_matrix_complex_fwrite
72 gsl_matrix_complex_fscanf
73 gsl_matrix_complex_fprintf
74 gsl_matrix_complex_memcpy
75 gsl_matrix_complex_swap
76 gsl_matrix_complex_swap_rows
77 gsl_matrix_complex_swap_columns
78 gsl_matrix_complex_swap_rowcol
79 gsl_matrix_complex_transpose
80 gsl_matrix_complex_transpose_memcpy
81 gsl_matrix_complex_isnull
82 gsl_matrix_complex_ispos
83 gsl_matrix_complex_isneg
84 gsl_matrix_complex_add
85 gsl_matrix_complex_sub
86 gsl_matrix_complex_mul_elements
87 gsl_matrix_complex_div_elements
88 gsl_matrix_complex_scale
89 gsl_matrix_complex_add_constant
90 gsl_matrix_complex_add_diagonal
91 gsl_matrix_complex_get_row
92 gsl_matrix_complex_get_col
93 gsl_matrix_complex_set_row
94 gsl_matrix_complex_set_col
98 my ($left,$right) = @_;
99 my $lcopy = $left->copy;
101 if ( blessed $right && $right->isa('Math::GSL::MatrixComplex') ) {
102 if ( $left->cols == $right->cols && $left->rows == $right->rows ) {
103 return _dot_product($left,$right);
104 #gsl_matrix_complex_mul_elements($lcopy->raw, $right->raw);
106 croak "Math::GSL::MatrixComplex - multiplication of elements of matrices must be called with two objects matrices and must have the same number of columns and rows";
109 gsl_matrix_complex_scale($lcopy->raw, $right);
115 my ($left,$right) = @_;
117 croak "Math::GSL::Matrix - incompatible matrices in multiplication"
118 unless ( blessed $right && $right->isa('Math::GSL::MatrixComplex') and $left->rows == $right->cols );
119 my $C = Math::GSL::MatrixComplex->new($left->rows, $right->cols);
120 my $complex = gsl_complex_rect(1,0);
121 gsl_blas_zgemm($CblasNoTrans, $CblasNoTrans, $complex, $left->raw, $right->raw, $complex, $C->raw);
126 my ($left, $right) = @_;
128 my $lcopy = $left->copy;
130 if ( blessed $right && $right->isa('Math::GSL::MatrixComplex') && blessed $left && $left->isa('Math::GSL::MatrixComplex') ) {
131 if ( $left->cols == $right->cols && $left->rows == $right->rows ) {
132 gsl_matrix_complex_add($lcopy->raw, $right->raw);
134 croak "Math::GSL - addition of matrices must be called with two objects matrices and must have the same number of columns and rows";
137 gsl_matrix_complex_add_constant($lcopy->raw, $right);
143 my ($left, $right) = @_;
145 my $lcopy = $left->copy;
147 if ( blessed $right && $right->isa('Math::GSL::MatrixComplex') && blessed $left && $left->isa('Math::GSL::MatrixComplex') ) {
148 if ( $left->cols == $right->cols && $left->rows == $right->rows ) {
149 gsl_matrix_complex_sub($lcopy->raw, $right->raw);
151 croak "Math::GSL - subtraction of matrices must be called with two objects matrices and must have the same number of columns and rows";
154 gsl_matrix_complex_add_constant($lcopy->raw, $right*-1);
160 my ($left, $right) = @_;
162 my $comp = Array::Compare->new;
163 return ($comp->compare([$left->as_list], [$right->as_list]));
167 my ($left, $right) = @_;
169 return !(&_equal($left, $right));
174 my $copy = Math::GSL::MatrixComplex->new( $self->rows, $self->cols );
175 if ( gsl_matrix_complex_memcpy($copy->raw, $self->raw) != $GSL_SUCCESS ) {
176 croak "Math::GSL - error copying memory, aborting";
180 %EXPORT_TAGS = ( all => \@EXPORT_OK );
184 Math::GSL::MatrixComplex - Complex Matrices
188 use Math::GSL::MatrixComplex qw/:all/;
189 my $matrix1 = Math::GSL::MatrixComplex->new(5,5); # OO interface
190 my $matrix3 = $matrix1 + $matrix1;
191 my $matrix4 = $matrix1 - $matrix1;
192 if($matrix1 == $matrix4) ...
193 if($matrix1 != $matrix3) ...
195 my $matrix2 = gsl_matrix_complex_alloc(5,5); # standard interface
197 =head1 Objected Oriented Interface to GSL Math::GSL::MatrixComplex
201 Creates a new MatrixComplex object of the given size.
203 my $matrix = Math::GSL::MatrixComplex->new(10,10);
208 my ($class, $rows, $cols) = @_;
211 if ( defined $rows && defined $cols &&
212 $rows > 0 && $cols > 0 &&
213 (int $rows == $rows) && (int $cols == $cols)){
215 $matrix = gsl_matrix_complex_alloc($rows,$cols);
217 croak( __PACKAGE__.'::new($x,$y) - $x and $y must be positive integers');
219 gsl_matrix_complex_set_zero($matrix);
220 $this->{_matrix} = $matrix;
221 ($this->{_rows}, $this->{_cols}) = ($rows,$cols);
226 Get the underlying GSL matrix object created by SWIG, useful for using gsl_matrix_* functions which do not have an OO counterpart.
228 my $matrix = Math::GSL::MatrixComplex->new(3,3);
229 my $gsl_matrix = $matrix->raw;
230 my $stuff = gsl_matrix_complex_get($gsl_matrix, 1, 2);
233 sub raw { (shift)->{_matrix} }
236 Returns the number of rows in the matrix.
238 my $rows = $matrix->rows;
241 sub rows { (shift)->{_rows} }
245 Returns the number of columns in the matrix.
247 my $cols = $matrix->cols;
250 sub cols { (shift)->{_cols} }
254 Returns a 1xN or Nx1 matrix as a Math::GSL::VectorComplex 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:
256 my $vector1 = $matrix->col(0)->as_vector;
257 my $vector2 = $matrix->row(1)->as_vector;
264 croak(__PACKAGE__.'::as_vector() - must be a single row or column matrix') unless ($self->cols == 1 or $self->rows == 1);
266 # TODO: there is a faster way to do this
267 return Math::GSL::VectorComplex->new([ $self->as_list ] );
273 Get the contents of a Math::GSL::Matrix object as a Perl list.
275 my $matrix = Math::GSL::MatrixComplex->new(3,3);
277 my @matrix = $matrix->as_list;
284 for my $row ( 0 .. $self->rows-1) {
285 push @matrix, map { gsl_matrix_complex_get($self->raw, $row, $_) } (0 .. $self->cols-1 );
296 Returns a row matrix of the row you enter.
298 my $matrix = Math::GSL::MatrixComplex->new(3,3);
300 my $matrix_row = $matrix->row(0);
306 my ($self, $row) = @_;
307 croak (__PACKAGE__.'::$matrix->row($row) - invalid $row value')
308 unless (($row < $self->rows) and $row >= 0);
310 my $rowmat = Math::GSL::MatrixComplex->new(1,$self->cols);
312 for my $n (0 .. $self->cols-1) {
313 gsl_matrix_complex_set( $rowmat->raw, 0, $n,
314 gsl_matrix_complex_get($self->raw, $row, $n)
323 Returns a col matrix of the column you enter.
325 my $matrix = Math::GSL::MatrixComplex->new(3,3);
327 my $matrix_col = $matrix->col(0);
333 my ($self, $col) = @_;
334 croak (__PACKAGE__."::\$matrix->col(\$col) - $col not a valid column")
335 unless (defined $col && $col < $self->cols and $col >= 0);
337 my $colmat = Math::GSL::MatrixComplex->new($self->rows, 1);
339 map { gsl_matrix_complex_set($colmat->raw, $_, 0,
340 gsl_matrix_complex_get($self->raw, $_, $col),
342 } (0 .. $self->rows-1);
349 Sets a the values of a row with the elements of an array.
351 my $matrix = Math::GSL::MatrixComplex->new(3,3);
352 $matrix->set_row(0, [8, 6, 2]);
354 You can also set multiple rows at once with chained calls:
355 my $matrix = Math::GSL::MatrixComplex->new(3,3);
356 $matrix->set_row(0, [8, 6, 2])
357 ->set_row(1, [2, 4, 1]);
363 my ($self, $row, $values) = @_;
364 my $length = $#$values;
365 die __PACKAGE__.'::set_row($x, $values) - $values must be a nonempty array reference' if $length == -1;
366 die __PACKAGE__.'::set_row($x, $values) - $x must be a valid row number' if ($row < 0 || $row >= $self->rows);
367 die __PACKAGE__.'::set_row($x, $values) - $values must contains the same number of elements as there is columns in the matrix' if($length != $self->cols-1);
368 # $values may have Math::Complex objects
369 @$values = map { Math::GSL::Complex::gsl_complex_rect(Re($_), Im($_)) } @$values;
370 # warn Dumper [ @$values ];
371 map { gsl_matrix_complex_set($self->raw, $row, $_, $values->[$_]) } (0..$length);
377 Sets a the values of a column with the elements of an array.
379 my $matrix = Math::GSL::MatrixComplex->new(3,3);
380 $matrix->set_col(0, [8, 6, 2]);
382 You can also set multiple columns at once with chained calls:
383 my $matrix = Math::GSL::MatrixComplex->new(3,3);
384 $matrix->set_col(0, [8, 6, 2])
385 ->set_col(1, [2, 4, 1]);
391 my ($self, $col, $values) = @_;
392 my $length = $#$values;
393 die __PACKAGE__.'::set_col($x, $values) - $values must be a nonempty array reference' if $length == -1;
394 die __PACKAGE__.'::set_col($x, $values) - $x must be a valid column number' if ($col < 0 || $col >= $self->cols);
395 die __PACKAGE__.'::set_col($x, $values) - $values must contains the same number of elements as there is rowss in the matrix' if($length != $self->rows-1);
396 # $values may have Math::Complex objects
397 @$values = map { gsl_complex_rect(Re($_), Im($_)) } @$values;
398 map { gsl_matrix_complex_set($self->raw, $_, $col, $values->[$_]) } (0..$length);
404 Returns true if a matrix is square, i.e. it has the same number of rows as columns, false otherwise.
411 return ($self->rows == $self->cols) ? 1 : 0 ;
416 Returns the determinant of a matrix (computed by LU decomposition) or dies if called on a non-square matrix.
418 my $det = $matrix->det();
425 croak(__PACKAGE__."- determinant only exists for square matrices") unless $self->is_square;
427 my $p = Math::GSL::Permutation->new( $self->rows );
428 my $LU = $self->copy;
430 my $s = gsl_linalg_complex_LU_decomp($LU->raw, $p->raw);
432 # $z is a gsl_complex
433 my $z = gsl_linalg_complex_LU_det($LU->raw, $s );
434 return Math::Complex->make( gsl_real($z), gsl_imag($z) );
439 Set a matrix to the zero matrix.
448 gsl_matrix_complex_set_zero($self->raw);
454 Set a matrix to the identity matrix, i.e. one on the diagonal and zero elsewhere.
456 my $I = $matrix->identity;
463 gsl_matrix_complex_set_identity($self->raw);
469 Returns the inverse of a matrix or dies when called on a non-square matrix.
471 my $inverse = $matrix->inverse;
478 croak(__PACKAGE__."- inverse only exists for square matrices") unless $self->is_square;
480 my $p = Math::GSL::Permutation->new( $self->rows );
481 my $LU = $self->copy;
482 my $inverse = $self->copy;
484 # should check return status
485 gsl_linalg_complex_LU_decomp($LU->raw, $p->raw);
486 gsl_linalg_complex_LU_invert($LU->raw, $p->raw,$inverse->raw);
492 =head2 is_hermitian()
494 Returns true if the matrix is hermitian, false otherwise
496 my $test = $matrix->is_hermitian;
503 my $test = $self->is_square;
504 my $transpose = $self->copy;
505 gsl_matrix_complex_transpose($transpose->raw);
508 for my $row (0..$self->rows - 1) {
509 map { gsl_matrix_complex_set($transpose->raw, $row, $_, gsl_complex_conjugate(gsl_matrix_complex_get($transpose->raw, $row, $_))) } (0..$self->cols - 1);
511 if($self != $transpose){
525 my ($r,$c) = ($self->rows,$self->cols);
526 croak "Math::GSL::MatrixComplex : \$matrix->eigenvalues - \$matrix must be square" unless ($r == $c);
527 my $evec = gsl_matrix_complex_alloc($r,$c);
528 my $eigen = gsl_eigen_nonsymmv_alloc($r);
529 my $vector = gsl_vector_complex_alloc($r);
532 # GSL has no generalized complex matrix routines
533 # can only continue if the matrix is hermitian
534 croak (__PACKAGE__."::eigenvalues : non-hermitian matrices are not currently supported") unless $self->is_hermitian;
535 gsl_eigen_nonsymmv($self->raw, $vector, $evec, $eigen);
537 my $x = gsl_vector_complex_real($vector); # vector of real components
538 my $y = gsl_vector_complex_imag($vector); # vector of imaginary components
541 _as_complex( gsl_vector_get($x->{vector}, $_) ,
542 gsl_vector_get($y->{vector}, $_) )
549 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.
551 my $lndet = $matrix->lndet();
558 croak(__PACKAGE__."- log determinant only exists for square matrices") unless $self->is_square;
560 my $p = Math::GSL::Permutation->new( $self->rows );
561 my $LU = $self->copy;
563 gsl_linalg_complex_LU_decomp($LU->raw, $p->raw);
564 return gsl_linalg_complex_LU_lndet($LU->raw);
570 =item C<gsl_matrix_complex_alloc >
572 =item C<gsl_matrix_complex_calloc >
574 =item C<gsl_matrix_complex_alloc_from_block >
576 =item C<gsl_matrix_complex_alloc_from_matrix >
578 =item C<gsl_vector_complex_alloc_row_from_matrix >
580 =item C<gsl_vector_complex_alloc_col_from_matrix >
582 =item C<gsl_matrix_complex_free >
584 =item C<gsl_matrix_complex_submatrix >
586 =item C<gsl_matrix_complex_row >
588 =item C<gsl_matrix_complex_column >
590 =item C<gsl_matrix_complex_diagonal >
592 =item C<gsl_matrix_complex_subdiagonal >
594 =item C<gsl_matrix_complex_superdiagonal >
596 =item C<gsl_matrix_complex_subrow >
598 =item C<gsl_matrix_complex_subcolumn >
600 =item C<gsl_matrix_complex_view_array >
602 =item C<gsl_matrix_complex_view_array_with_tda >
604 =item C<gsl_matrix_complex_view_vector >
606 =item C<gsl_matrix_complex_view_vector_with_tda >
608 =item C<gsl_matrix_complex_const_submatrix >
610 =item C<gsl_matrix_complex_const_row >
612 =item C<gsl_matrix_complex_const_column >
614 =item C<gsl_matrix_complex_const_diagonal >
616 =item C<gsl_matrix_complex_const_subdiagonal >
618 =item C<gsl_matrix_complex_const_superdiagonal >
620 =item C<gsl_matrix_complex_const_subrow >
622 =item C<gsl_matrix_complex_const_subcolumn >
624 =item C<gsl_matrix_complex_const_view_array >
626 =item C<gsl_matrix_complex_const_view_array_with_tda >
628 =item C<gsl_matrix_complex_const_view_vector >
630 =item C<gsl_matrix_complex_const_view_vector_with_tda >
632 =item C<gsl_matrix_complex_get>
634 =item C<gsl_matrix_complex_set>
636 =item C<gsl_matrix_complex_ptr>
638 =item C<gsl_matrix_complex_const_ptr>
640 =item C<gsl_matrix_complex_set_zero >
642 =item C<gsl_matrix_complex_set_identity >
644 =item C<gsl_matrix_complex_set_all >
646 =item C<gsl_matrix_complex_fread >
648 =item C<gsl_matrix_complex_fwrite >
650 =item C<gsl_matrix_complex_fscanf >
652 =item C<gsl_matrix_complex_fprintf >
654 =item C<gsl_matrix_complex_memcpy>
656 =item C<gsl_matrix_complex_swap>
658 =item C<gsl_matrix_complex_swap_rows>
660 =item C<gsl_matrix_complex_swap_columns>
662 =item C<gsl_matrix_complex_swap_rowcol>
664 =item C<gsl_matrix_complex_transpose >
666 =item C<gsl_matrix_complex_transpose_memcpy >
668 =item C<gsl_matrix_complex_isnull >
670 =item C<gsl_matrix_complex_ispos >
672 =item C<gsl_matrix_complex_isneg >
674 =item C<gsl_matrix_complex_add >
676 =item C<gsl_matrix_complex_sub >
678 =item C<gsl_matrix_complex_mul_elements >
680 =item C<gsl_matrix_complex_div_elements >
682 =item C<gsl_matrix_complex_scale >
684 =item C<gsl_matrix_complex_add_constant >
686 =item C<gsl_matrix_complex_add_diagonal >
688 =item C<gsl_matrix_complex_get_row>
690 =item C<gsl_matrix_complex_get_col>
692 =item C<gsl_matrix_complex_set_row>
694 =item C<gsl_matrix_complex_set_col>
698 For more informations on the functions, we refer you to the GSL offcial documentation
699 L<http://www.gnu.org/software/gsl/manual/html_node/>
704 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
706 =head1 COPYRIGHT AND LICENSE
708 Copyright (C) 2008 Jonathan Leto and Thierry Moisan
710 This program is free software; you can redistribute it and/or modify it
711 under the same terms as Perl itself.