3 use Scalar::Util 'blessed';
5 use Math::GSL qw/:all/;
6 use Math::GSL::Errno qw/:all/;
13 gsl_matrix_alloc gsl_matrix_calloc gsl_matrix_alloc_from_block
14 gsl_matrix_alloc_from_matrix gsl_vector_alloc_row_from_matrix
15 gsl_vector_alloc_col_from_matrix gsl_matrix_free gsl_matrix_submatrix
16 gsl_matrix_row gsl_matrix_column gsl_matrix_diagonal
17 gsl_matrix_subdiagonal gsl_matrix_superdiagonal gsl_matrix_subrow
18 gsl_matrix_subcolumn gsl_matrix_view_array
19 gsl_matrix_view_array_with_tda gsl_matrix_view_vector
20 gsl_matrix_view_vector_with_tda gsl_matrix_const_submatrix
21 gsl_matrix_const_row gsl_matrix_const_column gsl_matrix_const_diagonal
22 gsl_matrix_const_subdiagonal gsl_matrix_const_superdiagonal
23 gsl_matrix_const_subrow gsl_matrix_const_subcolumn
24 gsl_matrix_const_view_array gsl_matrix_const_view_array_with_tda
25 gsl_matrix_const_view_vector gsl_matrix_const_view_vector_with_tda
26 gsl_matrix_get gsl_matrix_set gsl_matrix_ptr gsl_matrix_const_ptr
27 gsl_matrix_set_zero gsl_matrix_set_identity gsl_matrix_set_all
28 gsl_matrix_fread gsl_matrix_fwrite gsl_matrix_fscanf gsl_matrix_fprintf
29 gsl_matrix_memcpy gsl_matrix_swap gsl_matrix_swap_rows
30 gsl_matrix_swap_columns gsl_matrix_swap_rowcol gsl_matrix_transpose
31 gsl_matrix_transpose_memcpy gsl_matrix_max gsl_matrix_min
32 gsl_matrix_minmax gsl_matrix_max_index gsl_matrix_min_index
33 gsl_matrix_minmax_index gsl_matrix_isnull gsl_matrix_ispos
34 gsl_matrix_isneg gsl_matrix_isnonneg gsl_matrix_add gsl_matrix_sub
35 gsl_matrix_mul_elements gsl_matrix_div_elements gsl_matrix_scale
36 gsl_matrix_add_constant gsl_matrix_add_diagonal
37 gsl_matrix_char_alloc gsl_matrix_char_calloc gsl_matrix_char_alloc_from_block
38 gsl_matrix_char_alloc_from_matrix gsl_vector_char_alloc_row_from_matrix gsl_vector_char_alloc_col_from_matrix
39 gsl_matrix_char_free gsl_matrix_char_submatrix
40 gsl_matrix_char_row gsl_matrix_char_column
41 gsl_matrix_char_diagonal gsl_matrix_char_subdiagonal gsl_matrix_char_superdiagonal
42 gsl_matrix_char_subrow gsl_matrix_char_subcolumn gsl_matrix_char_view_array
43 gsl_matrix_char_view_array_with_tda gsl_matrix_char_view_vector gsl_matrix_char_view_vector_with_tda
44 gsl_matrix_char_const_submatrix gsl_matrix_char_const_row gsl_matrix_char_const_column
45 gsl_matrix_char_const_diagonal gsl_matrix_char_const_subdiagonal gsl_matrix_char_const_superdiagonal
46 gsl_matrix_char_const_subrow gsl_matrix_char_const_subcolumn gsl_matrix_char_const_view_array
47 gsl_matrix_char_const_view_array_with_tda gsl_matrix_char_const_view_vector gsl_matrix_char_const_view_vector_with_tda
48 gsl_matrix_char_get gsl_matrix_char_set gsl_matrix_char_ptr gsl_matrix_char_const_ptr
49 gsl_matrix_char_set_zero gsl_matrix_char_set_identity
50 gsl_matrix_char_set_all gsl_matrix_char_fread
51 gsl_matrix_char_fwrite gsl_matrix_char_fscanf gsl_matrix_char_fprintf
52 gsl_matrix_char_memcpy gsl_matrix_char_swap
53 gsl_matrix_char_swap_rows gsl_matrix_char_swap_columns
54 gsl_matrix_char_swap_rowcol gsl_matrix_char_transpose gsl_matrix_char_transpose_memcpy
55 gsl_matrix_char_max gsl_matrix_char_min
56 gsl_matrix_char_minmax gsl_matrix_char_max_index
57 gsl_matrix_char_min_index gsl_matrix_char_minmax_index
58 gsl_matrix_char_isnull gsl_matrix_char_ispos gsl_matrix_char_isneg
59 gsl_matrix_char_isnonneg gsl_matrix_char_add
60 gsl_matrix_char_sub gsl_matrix_char_mul_elements gsl_matrix_char_div_elements
61 gsl_matrix_char_scale gsl_matrix_char_add_constant gsl_matrix_char_add_diagonal
62 gsl_matrix_int_alloc gsl_matrix_int_calloc gsl_matrix_int_alloc_from_block
63 gsl_matrix_int_alloc_from_matrix gsl_vector_int_alloc_row_from_matrix gsl_vector_int_alloc_col_from_matrix
64 gsl_matrix_int_free gsl_matrix_int_submatrix gsl_matrix_int_row
65 gsl_matrix_int_column gsl_matrix_int_diagonal gsl_matrix_int_subdiagonal
66 gsl_matrix_int_superdiagonal gsl_matrix_int_subrow gsl_matrix_int_subcolumn gsl_matrix_int_view_array
67 gsl_matrix_int_view_array_with_tda gsl_matrix_int_view_vector gsl_matrix_int_view_vector_with_tda
68 gsl_matrix_int_const_submatrix gsl_matrix_int_const_row gsl_matrix_int_const_column
69 gsl_matrix_int_const_diagonal gsl_matrix_int_const_subdiagonal gsl_matrix_int_const_superdiagonal
70 gsl_matrix_int_const_subrow gsl_matrix_int_const_subcolumn gsl_matrix_int_const_view_array
71 gsl_matrix_int_const_view_array_with_tda gsl_matrix_int_const_view_vector gsl_matrix_int_const_view_vector_with_tda
72 gsl_matrix_int_get gsl_matrix_int_set
73 gsl_matrix_int_ptr gsl_matrix_int_const_ptr
74 gsl_matrix_int_set_zero gsl_matrix_int_set_identity gsl_matrix_int_set_all
75 gsl_matrix_int_fread gsl_matrix_int_fwrite
76 gsl_matrix_int_fscanf gsl_matrix_int_fprintf
77 gsl_matrix_int_memcpy gsl_matrix_int_swap
78 gsl_matrix_int_swap_rows gsl_matrix_int_swap_columns gsl_matrix_int_swap_rowcol
79 gsl_matrix_int_transpose gsl_matrix_int_transpose_memcpy
80 gsl_matrix_int_max gsl_matrix_int_min gsl_matrix_int_minmax
81 gsl_matrix_int_max_index gsl_matrix_int_min_index
82 gsl_matrix_int_minmax_index gsl_matrix_int_isnull
83 gsl_matrix_int_ispos gsl_matrix_int_isneg gsl_matrix_int_isnonneg
84 gsl_matrix_int_add gsl_matrix_int_sub
85 gsl_matrix_int_mul_elements gsl_matrix_int_div_elements gsl_matrix_int_scale
86 gsl_matrix_int_add_constant gsl_matrix_int_add_diagonal
87 gsl_matrix_get_row gsl_matrix_get_col gsl_matrix_set_row gsl_matrix_set_col
88 gsl_matrix_complex_alloc gsl_matrix_complex_calloc gsl_matrix_complex_alloc_from_block
89 gsl_matrix_complex_alloc_from_matrix gsl_vector_complex_alloc_row_from_matrix gsl_vector_complex_alloc_col_from_matrix
90 gsl_matrix_complex_free gsl_matrix_complex_submatrix gsl_matrix_complex_row
91 gsl_matrix_complex_column gsl_matrix_complex_diagonal gsl_matrix_complex_subdiagonal
92 gsl_matrix_complex_superdiagonal gsl_matrix_complex_subrow gsl_matrix_complex_subcolumn
93 gsl_matrix_complex_view_array gsl_matrix_complex_view_array_with_tda gsl_matrix_complex_view_vector
94 gsl_matrix_complex_view_vector_with_tda gsl_matrix_complex_const_submatrix gsl_matrix_complex_const_row
95 gsl_matrix_complex_const_column gsl_matrix_complex_const_diagonal gsl_matrix_complex_const_subdiagonal
96 gsl_matrix_complex_const_superdiagonal gsl_matrix_complex_const_subrow gsl_matrix_complex_const_subcolumn
97 gsl_matrix_complex_const_view_array gsl_matrix_complex_const_view_array_with_tda gsl_matrix_complex_const_view_vector
98 gsl_matrix_complex_const_view_vector_with_tda gsl_matrix_complex_get gsl_matrix_complex_set
99 gsl_matrix_complex_ptr gsl_matrix_complex_const_ptr gsl_matrix_complex_set_zero
100 gsl_matrix_complex_set_identity gsl_matrix_complex_set_all gsl_matrix_complex_fread
101 gsl_matrix_complex_fwrite gsl_matrix_complex_fscanf gsl_matrix_complex_fprintf
102 gsl_matrix_complex_memcpy gsl_matrix_complex_swap gsl_matrix_complex_swap_rows
103 gsl_matrix_complex_swap_columns gsl_matrix_complex_swap_rowcol gsl_matrix_complex_transpose
104 gsl_matrix_complex_transpose_memcpy gsl_matrix_complex_isnull gsl_matrix_complex_ispos
105 gsl_matrix_complex_isneg gsl_matrix_complex_add gsl_matrix_complex_sub
106 gsl_matrix_complex_mul_elements gsl_matrix_complex_div_elements gsl_matrix_complex_scale
107 gsl_matrix_complex_add_constant gsl_matrix_complex_add_diagonal gsl_matrix_complex_get_row
108 gsl_matrix_complex_get_col gsl_matrix_complex_set_row gsl_matrix_complex_set_col /;
112 %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
276 gsl_matrix_complex_alloc
277 gsl_matrix_complex_calloc
278 gsl_matrix_complex_alloc_from_block
279 gsl_matrix_complex_alloc_from_matrix
280 gsl_vector_complex_alloc_row_from_matrix
281 gsl_vector_complex_alloc_col_from_matrix
282 gsl_matrix_complex_free
283 gsl_matrix_complex_submatrix
284 gsl_matrix_complex_row
285 gsl_matrix_complex_column
286 gsl_matrix_complex_diagonal
287 gsl_matrix_complex_subdiagonal
288 gsl_matrix_complex_superdiagonal
289 gsl_matrix_complex_subrow
290 gsl_matrix_complex_subcolumn
291 gsl_matrix_complex_view_array
292 gsl_matrix_complex_view_array_with_tda
293 gsl_matrix_complex_view_vector
294 gsl_matrix_complex_view_vector_with_tda
295 gsl_matrix_complex_const_submatrix
296 gsl_matrix_complex_const_row
297 gsl_matrix_complex_const_column
298 gsl_matrix_complex_const_diagonal
299 gsl_matrix_complex_const_subdiagonal
300 gsl_matrix_complex_const_superdiagonal
301 gsl_matrix_complex_const_subrow
302 gsl_matrix_complex_const_subcolumn
303 gsl_matrix_complex_const_view_array
304 gsl_matrix_complex_const_view_array_with_tda
305 gsl_matrix_complex_const_view_vector
306 gsl_matrix_complex_const_view_vector_with_tda
307 gsl_matrix_complex_get
308 gsl_matrix_complex_set
309 gsl_matrix_complex_ptr
310 gsl_matrix_complex_const_ptr
311 gsl_matrix_complex_set_zero gsl_matrix_complex_set_identity
312 gsl_matrix_complex_set_all
313 gsl_matrix_complex_fread
314 gsl_matrix_complex_fwrite
315 gsl_matrix_complex_fscanf
316 gsl_matrix_complex_fprintf
317 gsl_matrix_complex_memcpygsl_matrix_complex_swap
318 gsl_matrix_complex_swap_rows
319 gsl_matrix_complex_swap_columns
320 gsl_matrix_complex_swap_rowcol
321 gsl_matrix_complex_transpose
322 gsl_matrix_complex_transpose_memcpy
323 gsl_matrix_complex_isnull
324 gsl_matrix_complex_ispos
325 gsl_matrix_complex_isneg
326 gsl_matrix_complex_add
327 gsl_matrix_complex_sub
328 gsl_matrix_complex_mul_elements
329 gsl_matrix_complex_div_elements
330 gsl_matrix_complex_scale
331 gsl_matrix_complex_add_constant
332 gsl_matrix_complex_add_diagonal
333 gsl_matrix_complex_get_row
334 gsl_matrix_complex_get_col
335 gsl_matrix_complex_set_row
336 gsl_matrix_complex_set_col
343 Math::GSL::Matrix - Mathematical functions concerning Matrices
347 use Math::GSL::Matrix qw/:all/;
348 my $matrix1 = Math::GSL::Matrix->new(5,5); # OO interface
349 my $matrix2 = $matrix1 + 4; # You can add or substract values or matrices to OO matrices
350 my $matrix3 = $matrix1 - 4;
351 my $matrix4 = $matrix2 + $matrix1;
352 my $matrix5 = gsl_matrix_alloc(5,5); # standard interface
355 =head1 Objected Oriented Interface to GSL Math::GSL::Matrix
357 =head2 Math::GSL::Matrix->new()
359 Creates a new Matrix of the given size.
361 my $matrix = Math::GSL::Matrix->new(10,10);
366 my ($class, $rows, $cols) = @_;
369 if ( defined $rows && defined $cols &&
370 $rows > 0 && $cols > 0 &&
371 (int $rows == $rows) && (int $cols == $cols)){
373 $matrix = gsl_matrix_alloc($rows,$cols);
375 croak( __PACKAGE__.'::new($x,$y) - $x and $y must be positive integers');
377 gsl_matrix_set_zero($matrix);
378 $this->{_matrix} = $matrix;
379 ($this->{_rows}, $this->{_cols}) = ($rows,$cols);
384 Get the underlying GSL matrix object created by SWIG, useful for using gsl_matrix_* functions which do not have an OO counterpart.
386 my $matrix = Math::GSL::Matrix->new(3,3);
387 my $gsl_matrix = $matrix->raw;
388 my $stuff = gsl_matrix_get($gsl_matrix, 1, 2);
391 sub raw { (shift)->{_matrix} }
395 Returns a copy of the matrix, which has the same size and values but resides at a different location in memory.
397 my $matrix = Math::GSL::Matrix->new(5,5);
398 my $copy = $matrix->copy;
404 my $copy = Math::GSL::Matrix->new( $self->rows, $self->cols );
405 if ( gsl_matrix_memcpy($copy->raw, $self->raw) != $GSL_SUCCESS ) {
406 croak "Math::GSL - error copying memory, aborting";
414 Returns the number of rows in the matrix.
416 my $rows = $matrix->rows;
419 sub rows { (shift)->{_rows} }
423 Returns the number of columns in the matrix.
425 my $cols = $matrix->cols;
428 sub cols { (shift)->{_cols} }
432 Get the contents of a Math::GSL::Matrix object as a Perl list.
434 my $matrix = Math::GSL::Matrix->new(3,3);
436 my @matrix = $matrix->as_list;
446 for($line=0; $line<$self->rows; $line++){
448 gsl_matrix_get($self->raw, $line, $_)
449 } (0 .. $self->cols-1 );
457 Returns a row matrix of the row you enter.
459 my $matrix = Math::GSL::Matrix->new(3,3);
461 my $matrix_row = $matrix->row(0);
467 my ($self, $row) = @_;
468 croak (__PACKAGE__.'::$matrix->row($row) - invalid $row value')
469 unless (($row < $self->rows) and $row >= 0);
471 my $rowmat = Math::GSL::Matrix->new(1,$self->cols);
474 for my $n (0 .. $self->cols-1) {
475 push (@got, gsl_matrix_get($self->raw, $row, $n));
477 for my $n (0 .. $self->cols-1) {
478 gsl_matrix_set($rowmat->raw, 0, $n, $got[$n]); }
485 Returns a col matrix of the column you enter.
487 my $matrix = Math::GSL::Matrix->new(3,3);
489 my $matrix_col = $matrix->col(0);
495 my ($self, $col) = @_;
496 croak (__PACKAGE__."::\$matrix->col(\$col) - $col not a valid column")
497 unless ($col < $self->cols and $col >= 0);
499 my $colvec = Math::GSL::Vector->new($self->cols);
500 my $colmat = Math::GSL::Matrix->new($self->rows, 1);
503 for my $n (0 .. $self->rows-1) {
504 push (@got, gsl_matrix_get($self->raw, $n, $col));
506 for my $n (0 .. $self->rows-1) {
507 gsl_matrix_set($colmat->raw, $n, 0, $got[$n]); }
514 Sets a the values of a row with the elements of an array.
516 my $matrix = Math::GSL::Matrix->new(3,3);
517 $matrix->set_row(0, [8, 6, 2]);
519 You can also set multiple rows at once with chained calls:
520 my $matrix = Math::GSL::Matrix->new(3,3);
521 $matrix->set_row(0, [8, 6, 2])
522 ->set_row(1, [2, 4, 1]);
528 my ($self, $row, $values) = @_;
529 my $length = $#$values;
530 die __PACKAGE__.'::new($x, $values) - $values must be a nonempty array reference' if $length == -1;
531 die __PACKAGE__.'::set_row($x, $values) - $x must be a valid row number' if ($row < 0 || $row >= $self->rows);
532 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);
533 map { gsl_matrix_set($self->raw, $row, $_, $values->[$_]) } (0..$length);
539 Sets a the values of a column with the elements of an array.
541 my $matrix = Math::GSL::Matrix->new(3,3);
542 $matrix->set_col(0, [8, 6, 2]);
544 You can also set multiple columns at once with chained calls:
545 my $matrix = Math::GSL::Matrix->new(3,3);
546 $matrix->set_col(0, [8, 6, 2])
547 ->set_col(1, [2, 4, 1]);
553 my ($self, $col, $values) = @_;
554 my $length = $#$values;
555 die __PACKAGE__.'::new($x, $values) - $values must be a nonempty array reference' if $length == -1;
556 die __PACKAGE__.'::set_row($x, $values) - $x must be a valid column number' if ($col < 0 || $col >= $self->cols);
557 die __PACKAGE__.'::set_row($x, $values) - $values must contains the same number of elements as there is rowss in the matrix' if($length != $self->rows-1);
558 map { gsl_matrix_set($self->raw, $_, $col, $values->[$_]) } (0..$length);
563 my ($left, $right) = @_;
565 my $lcopy = $left->copy;
567 if ( blessed $right && $right->isa('Math::GSL::Matrix') && blessed $left && $left->isa('Math::GSL::Matrix') ) {
568 if ( $left->cols == $right->cols && $left->rows == $right->rows ) {
569 gsl_matrix_add($lcopy->raw, $right->raw);
571 croak "Math::GSL - addition of matrices must be called with two objects matrices and must have the same number of columnss and rows";
574 gsl_matrix_add_constant($lcopy->raw, $right);
580 my ($left, $right) = @_;
582 my $lcopy = $left->copy;
584 if ( blessed $right && $right->isa('Math::GSL::Matrix') && blessed $left && $left->isa('Math::GSL::Matrix') ) {
585 if ( $left->cols == $right->cols && $left->rows == $right->rows ) {
586 gsl_matrix_sub($lcopy->raw, $right->raw);
588 croak "Math::GSL - subtraction of matrices must be called with two objects matrices and must have the same number of columnss and rows";
591 gsl_matrix_add_constant($lcopy->raw, $right*-1);
598 Here is a list of all the functions included in this module :
602 =item C<gsl_matrix_alloc($i, $j)> - Return a gsl_matrix of $i rows and $j columns
604 =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
606 =item C<gsl_matrix_alloc_from_block> -
608 =item C<gsl_matrix_free> -
610 =item C<gsl_matrix_alloc_from_matrix > -
612 =item C<gsl_vector_alloc_row_from_matrix> -
614 =item C<gsl_vector_alloc_col_from_matrix > -
616 =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.
618 =item C<gsl_matrix_row($m , $i)> - Return a vector view of the $i-th row of the matrix $m
620 =item C<gsl_matrix_column($m, $j)> - Return a vector view of the $j-th column of the matrix $m
622 =item C<gsl_matrix_diagonal($m)> - Return a vector view of the diagonal of the vector. The matrix doesn't have to be square.
624 =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.
626 =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.
628 =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.
630 =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.
632 =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.
634 =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.
636 =item C<gsl_matrix_view_vector> -
638 =item C<gsl_matrix_view_vector_with_tda> -
640 =item C<gsl_matrix_const_submatrix> -
642 =item C<gsl_matrix_get($m, $i, $j)> - Return the (i,j)-th element of the matrix $m
644 =item C<gsl_matrix_set($m, $i, $j, $x)> - Set the value of the (i,j)-th element of the matrix $m to $x
646 =item C<gsl_matrix_ptr> -
648 =item C<gsl_matrix_const_ptr> -
650 =item C<gsl_matrix_set_zero($m)> - Set all the elements of the matrix $m to zero
652 =item C<gsl_matrix_set_identity($m)> - Set the elements of the matrix $m to the corresponding elements of the identity matrix
654 =item C<gsl_matrix_set_all($m, $x)> - Set all the elements of the matrix $m to the value $x
656 =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
658 =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
660 =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
662 =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
664 =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.
666 =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.
668 =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.
670 =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.
672 =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.
674 =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.
676 =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.
678 =item C<gsl_matrix_max($m)> - Return the maximum value in the matrix $m
680 =item C<gsl_matrix_min($m)> - Return the minimum value in the matrix $m
682 =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.
684 =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.
686 =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.
688 =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
690 =item C<gsl_matrix_isnull($m)> - Return 1 if all the elements of the matrix $m are zero, 0 otherwise
692 =item C<gsl_matrix_ispos($m)> - Return 1 if all the elements of the matrix $m are strictly positve, 0 otherwise
694 =item C<gsl_matrix_isneg($m)> - Return 1 if all the elements of the matrix $m are strictly negative, 0 otherwise
696 =item C<gsl_matrix_isnonneg($m)> - Return 1 if all the elements of the matrix $m are non-negatuive, 0 otherwise
698 =item C<gsl_matrix_add($a, $b)> - Add the elements of matrix $b to the elements of matrix $a
700 =item C<gsl_matrix_sub($a, $b)> - Subtract the elements of matrix $b from the elements of matrix $a
702 =item C<gsl_matrix_mul_elements($a, $b)> - Multiplie the elements of matrix $a by the elements of matrix $b
704 =item C<gsl_matrix_div_elements($a, $b)> - Divide the elements of matrix $a by the elements of matrix $b
706 =item C<gsl_matrix_scale($a, $x)> - Multiplie the elements of matrix $a by the constant factor $x
708 =item C<gsl_matrix_add_constant($a, $x)> - Add the constant value $x to the elements of the matrix $a
710 =item C<gsl_matrix_add_diagonal($a, $x)> - Add the constant value $x to the elements of the diagonal of the matrix $a
712 =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.
714 =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.
716 =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.
718 =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.
722 The following functions are specific to matrices containing complex numbers :
726 =item C<gsl_matrix_complex_alloc >
728 =item C<gsl_matrix_complex_calloc >
730 =item C<gsl_matrix_complex_alloc_from_block >
732 =item C<gsl_matrix_complex_alloc_from_matrix >
734 =item C<gsl_vector_complex_alloc_row_from_matrix >
736 =item C<gsl_vector_complex_alloc_col_from_matrix >
738 =item C<gsl_matrix_complex_free >
740 =item C<gsl_matrix_complex_submatrix >
742 =item C<gsl_matrix_complex_row >
744 =item C<gsl_matrix_complex_column >
746 =item C<gsl_matrix_complex_diagonal >
748 =item C<gsl_matrix_complex_subdiagonal >
750 =item C<gsl_matrix_complex_superdiagonal >
752 =item C<gsl_matrix_complex_subrow >
754 =item C<gsl_matrix_complex_subcolumn >
756 =item C<gsl_matrix_complex_view_array >
758 =item C<gsl_matrix_complex_view_array_with_tda >
760 =item C<gsl_matrix_complex_view_vector >
762 =item C<gsl_matrix_complex_view_vector_with_tda >
764 =item C<gsl_matrix_complex_const_submatrix >
766 =item C<gsl_matrix_complex_const_row >
768 =item C<gsl_matrix_complex_const_column >
770 =item C<gsl_matrix_complex_const_diagonal >
772 =item C<gsl_matrix_complex_const_subdiagonal >
774 =item C<gsl_matrix_complex_const_superdiagonal >
776 =item C<gsl_matrix_complex_const_subrow >
778 =item C<gsl_matrix_complex_const_subcolumn >
780 =item C<gsl_matrix_complex_const_view_array >
782 =item C<gsl_matrix_complex_const_view_array_with_tda >
784 =item C<gsl_matrix_complex_const_view_vector >
786 =item C<gsl_matrix_complex_const_view_vector_with_tda >
788 =item C<gsl_matrix_complex_get>
790 =item C<gsl_matrix_complex_set>
792 =item C<gsl_matrix_complex_ptr>
794 =item C<gsl_matrix_complex_const_ptr>
796 =item C<gsl_matrix_complex_set_zero >
798 =item C<gsl_matrix_complex_set_identity >
800 =item C<gsl_matrix_complex_set_all >
802 =item C<gsl_matrix_complex_fread >
804 =item C<gsl_matrix_complex_fwrite >
806 =item C<gsl_matrix_complex_fscanf >
808 =item C<gsl_matrix_complex_fprintf >
810 =item C<gsl_matrix_complex_memcpy>
812 =item C<gsl_matrix_complex_swap>
814 =item C<gsl_matrix_complex_swap_rows>
816 =item C<gsl_matrix_complex_swap_columns>
818 =item C<gsl_matrix_complex_swap_rowcol>
820 =item C<gsl_matrix_complex_transpose >
822 =item C<gsl_matrix_complex_transpose_memcpy >
824 =item C<gsl_matrix_complex_isnull >
826 =item C<gsl_matrix_complex_ispos >
828 =item C<gsl_matrix_complex_isneg >
830 =item C<gsl_matrix_complex_add >
832 =item C<gsl_matrix_complex_sub >
834 =item C<gsl_matrix_complex_mul_elements >
836 =item C<gsl_matrix_complex_div_elements >
838 =item C<gsl_matrix_complex_scale >
840 =item C<gsl_matrix_complex_add_constant >
842 =item C<gsl_matrix_complex_add_diagonal >
844 =item C<gsl_matrix_complex_get_row>
846 =item C<gsl_matrix_complex_get_col>
848 =item C<gsl_matrix_complex_set_row>
850 =item C<gsl_matrix_complex_set_col>
854 These are related to constant views of a matrix.
859 =item C<gsl_matrix_const_row>
861 =item C<gsl_matrix_const_colum>
863 =item C<gsl_matrix_const_diagonal>
865 =item C<gsl_matrix_const_subdiagonal>
867 =item C<gsl_matrix_const_superdiagonal>
869 =item C<gsl_matrix_const_subrow>
871 =item C<gsl_matrix_const_subcolumn>
873 =item C<gsl_matrix_const_view_array>
875 =item C<gsl_matrix_const_view_array_with_tda>
880 The following functions are similar to those above but work with C<char>'s and C<int>'s. We are not quite
881 sure if anyone wants these. Please speak up if you do and/or submit some patches to this documentation, please!
886 =item gsl_matrix_const_view_vector
888 =item gsl_matrix_const_view_vector_with_tda
890 =item gsl_matrix_char_alloc
892 =item gsl_matrix_char_calloc
894 =item gsl_matrix_char_alloc_from_block
896 =item gsl_matrix_char_alloc_from_matrix
898 =item gsl_vector_char_alloc_row_from_matrix
900 =item gsl_vector_char_alloc_col_from_matrix
902 =item gsl_matrix_char_free
904 =item gsl_matrix_char_submatrix
906 =item gsl_matrix_char_row
908 =item gsl_matrix_char_column
910 =item gsl_matrix_char_diagonal
912 =item gsl_matrix_char_subdiagonal
914 =item gsl_matrix_char_superdiagonal
916 =item gsl_matrix_char_subrow
918 =item gsl_matrix_char_subcolumn
920 =item gsl_matrix_char_view_array
922 =item gsl_matrix_char_view_array_with_tda
924 =item gsl_matrix_char_view_vector
926 =item gsl_matrix_char_view_vector_with_tda
928 =item gsl_matrix_char_const_submatrix
930 =item gsl_matrix_char_const_row
932 =item gsl_matrix_char_const_column
934 =item gsl_matrix_char_const_diagonal
936 =item gsl_matrix_char_const_subdiagonal
938 =item gsl_matrix_char_const_superdiagonal
940 =item gsl_matrix_char_const_subrow
942 =item gsl_matrix_char_const_subcolumn
944 =item gsl_matrix_char_const_view_array
946 =item gsl_matrix_char_const_view_array_with_tda
948 =item gsl_matrix_char_const_view_vector
950 =item gsl_matrix_char_const_view_vector_with_tda
952 =item gsl_matrix_char_get
954 =item gsl_matrix_char_set
956 =item gsl_matrix_char_ptr
958 =item gsl_matrix_char_const_ptr
960 =item gsl_matrix_char_set_zero
962 =item gsl_matrix_char_set_identity
964 =item gsl_matrix_char_set_all
966 =item gsl_matrix_char_fread
968 =item gsl_matrix_char_fwrite
970 =item gsl_matrix_char_fscanf
972 =item gsl_matrix_char_fprintf
974 =item gsl_matrix_char_memcpy
976 =item gsl_matrix_char_swap
978 =item gsl_matrix_char_swap_rows
980 =item gsl_matrix_char_swap_columns
982 =item gsl_matrix_char_swap_rowcol
984 =item gsl_matrix_char_transpose
986 =item gsl_matrix_char_transpose_memcpy
988 =item gsl_matrix_char_max
990 =item gsl_matrix_char_min
992 =item gsl_matrix_char_minmax
994 =item gsl_matrix_char_max_index
996 =item gsl_matrix_char_min_index
998 =item gsl_matrix_char_minmax_index
1000 =item gsl_matrix_char_isnull
1002 =item gsl_matrix_char_ispos
1004 =item gsl_matrix_char_isneg
1006 =item gsl_matrix_char_isnonneg
1008 =item gsl_matrix_char_add
1010 =item gsl_matrix_char_sub
1012 =item gsl_matrix_char_mul_elements
1014 =item gsl_matrix_char_div_elements
1016 =item gsl_matrix_char_scale
1018 =item gsl_matrix_char_add_constant
1020 =item gsl_matrix_char_add_diagonal
1022 =item gsl_matrix_int_alloc
1024 =item gsl_matrix_int_calloc
1026 =item gsl_matrix_int_alloc_from_block
1028 =item gsl_matrix_int_alloc_from_matrix
1030 =item gsl_vector_int_alloc_row_from_matrix
1032 =item gsl_vector_int_alloc_col_from_matrix
1034 =item gsl_matrix_int_free
1036 =item gsl_matrix_int_submatrix
1038 =item gsl_matrix_int_row
1040 =item gsl_matrix_int_column
1042 =item gsl_matrix_int_diagonal
1044 =item gsl_matrix_int_subdiagonal
1046 =item gsl_matrix_int_superdiagonal
1048 =item gsl_matrix_int_subrow
1050 =item gsl_matrix_int_subcolumn
1052 =item gsl_matrix_int_view_array
1054 =item gsl_matrix_int_view_array_with_tda
1056 =item gsl_matrix_int_view_vector
1058 =item gsl_matrix_int_view_vector_with_tda
1060 =item gsl_matrix_int_const_submatrix
1062 =item gsl_matrix_int_const_row
1064 =item gsl_matrix_int_const_column
1066 =item gsl_matrix_int_ptr
1068 =item gsl_matrix_int_const_ptr
1070 =item gsl_matrix_int_set_zero
1072 =item gsl_matrix_int_set_identity
1074 =item gsl_matrix_int_set_all
1076 =item gsl_matrix_int_fread
1078 =item gsl_matrix_int_fwrite
1080 =item gsl_matrix_int_fscanf
1082 =item gsl_matrix_int_fprintf
1084 =item gsl_matrix_int_memcpy
1086 =item gsl_matrix_int_swap
1088 =item gsl_matrix_int_swap_rows
1090 =item gsl_matrix_int_swap_columns
1092 =item gsl_matrix_int_swap_rowcol
1094 =item gsl_matrix_int_transpose
1096 =item gsl_matrix_int_transpose_memcpy
1098 =item gsl_matrix_int_max
1100 =item gsl_matrix_int_min
1102 =item gsl_matrix_int_minmax
1104 =item gsl_matrix_int_max_index
1106 =item gsl_matrix_int_min_index
1108 =item gsl_matrix_int_minmax_index
1110 =item gsl_matrix_int_isnull
1112 =item gsl_matrix_int_ispos
1114 =item gsl_matrix_int_isneg
1116 =item gsl_matrix_int_isnonneg
1118 =item gsl_matrix_int_add
1120 =item gsl_matrix_int_sub
1122 =item gsl_matrix_int_mul_elements
1124 =item gsl_matrix_int_div_elements
1126 =item gsl_matrix_int_scale
1128 =item gsl_matrix_int_add_constant
1130 =item gsl_matrix_int_add_diagonal
1134 You have to add the functions you want to use inside the qw /put_funtion_here /.
1135 You can also write use Math::GSL::Matrix qw/:all/ to use all avaible functions of the module.
1136 Other tags are also avaible, here is a complete list of all tags for this module :
1152 For more informations on the functions, we refer you to the GSL offcial documentation
1153 L<http://www.gnu.org/software/gsl/manual/html_node/>
1155 Tip : search on google: site:http://www.gnu.org/software/gsl/manual/html_node/name_of_the_function_you_want
1159 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>
1161 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.
1163 use Math::GSL::Matrix qw/:all/;
1164 my $m = gsl_matrix_alloc (10,3);
1167 gsl_matrix_set($m, $i, $j, 0.23 + 100*$i + $j);
1171 for my $i (0..99){ # OUT OF RANGE ERROR
1173 print "m($i, $j) = " . gsl_matrix_get ($m, $i, $j) . "\n";
1176 gsl_matrix_free ($m);
1179 use Math::GSL::Matrix qw/:all/;
1181 my $m = gsl_matrix_alloc (100, 100);
1182 my $a = gsl_matrix_alloc (100, 100);
1186 gsl_matrix_set ($m, $i, $j, 0.23 + $i + $j);
1190 The next program shows how to write a matrix to a file.
1192 my $out = gsl_fopen("test.dat", "wb");
1193 gsl_matrix_fwrite ($out, $m);
1196 my $in = gsl_fopen("test.dat", "rb");
1197 gsl_matrix_fread ($in, $a);
1203 $mij = gsl_matrix_get ($m, $i, $j);
1204 $aij = gsl_matrix_get ($a, $i, $j);
1205 $k++ if ($mij != $aij);
1209 gsl_matrix_free($m);
1210 gsl_matrix_free($a);
1212 print "differences = $k (should be zero)\n";
1216 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
1218 =head1 COPYRIGHT AND LICENSE
1220 Copyright (C) 2008 Jonathan Leto and Thierry Moisan
1222 This program is free software; you can redistribute it and/or modify it
1223 under the same terms as Perl itself.