Sort numerically, not asciibetically
[Math-GSL.git] / Matrix.i
blob049da1ca076e0a9425e78955e45a2b31517b847d
1 %module "Math::GSL::Matrix"
3 %include "typemaps.i"
4 %include "gsl_typemaps.i"
5 %apply int *OUTPUT { size_t *imin, size_t *imax, size_t *jmin, size_t *jmax };
6 %apply double *OUTPUT { double * min_out, double * max_out };
8 %{
9 #include "gsl/gsl_matrix.h"
10 #include "gsl/gsl_complex.h"
11 #include "gsl/gsl_vector_double.h"
12 #include "gsl/gsl_matrix_double.h"
13 #include "gsl/gsl_matrix_int.h"
14 #include "gsl/gsl_matrix_complex_double.h"
15 #include "gsl/gsl_matrix_char.h"
18 %include "gsl/gsl_matrix.h"
19 %include "gsl/gsl_complex.h"
20 %include "gsl/gsl_vector_double.h"
21 %include "gsl/gsl_matrix_double.h"
22 %include "gsl/gsl_matrix_int.h"
23 %include "gsl/gsl_matrix_complex_double.h"
24 %include "gsl/gsl_matrix_char.h"
27 %perlcode %{
29 use Carp qw/croak/;
30 use Math::GSL qw/:all/;
31 use Math::GSL::Errno qw/:all/;
33 @EXPORT_OK = qw/
34 gsl_matrix_alloc gsl_matrix_calloc gsl_matrix_alloc_from_block
35 gsl_matrix_alloc_from_matrix gsl_vector_alloc_row_from_matrix
36 gsl_vector_alloc_col_from_matrix gsl_matrix_free gsl_matrix_submatrix
37 gsl_matrix_row gsl_matrix_column gsl_matrix_diagonal
38 gsl_matrix_subdiagonal gsl_matrix_superdiagonal gsl_matrix_subrow
39 gsl_matrix_subcolumn gsl_matrix_view_array
40 gsl_matrix_view_array_with_tda gsl_matrix_view_vector
41 gsl_matrix_view_vector_with_tda gsl_matrix_const_submatrix
42 gsl_matrix_const_row gsl_matrix_const_column gsl_matrix_const_diagonal
43 gsl_matrix_const_subdiagonal gsl_matrix_const_superdiagonal
44 gsl_matrix_const_subrow gsl_matrix_const_subcolumn
45 gsl_matrix_const_view_array gsl_matrix_const_view_array_with_tda
46 gsl_matrix_const_view_vector gsl_matrix_const_view_vector_with_tda
47 gsl_matrix_get gsl_matrix_set gsl_matrix_ptr gsl_matrix_const_ptr
48 gsl_matrix_set_zero gsl_matrix_set_identity gsl_matrix_set_all
49 gsl_matrix_fread gsl_matrix_fwrite gsl_matrix_fscanf gsl_matrix_fprintf
50 gsl_matrix_memcpy gsl_matrix_swap gsl_matrix_swap_rows
51 gsl_matrix_swap_columns gsl_matrix_swap_rowcol gsl_matrix_transpose
52 gsl_matrix_transpose_memcpy gsl_matrix_max gsl_matrix_min
53 gsl_matrix_minmax gsl_matrix_max_index gsl_matrix_min_index
54 gsl_matrix_minmax_index gsl_matrix_isnull gsl_matrix_ispos
55 gsl_matrix_isneg gsl_matrix_isnonneg gsl_matrix_add gsl_matrix_sub
56 gsl_matrix_mul_elements gsl_matrix_div_elements gsl_matrix_scale
57 gsl_matrix_add_constant gsl_matrix_add_diagonal
58 gsl_matrix_char_alloc gsl_matrix_char_calloc gsl_matrix_char_alloc_from_block
59 gsl_matrix_char_alloc_from_matrix gsl_vector_char_alloc_row_from_matrix gsl_vector_char_alloc_col_from_matrix
60 gsl_matrix_char_free gsl_matrix_char_submatrix
61 gsl_matrix_char_row gsl_matrix_char_column
62 gsl_matrix_char_diagonal gsl_matrix_char_subdiagonal gsl_matrix_char_superdiagonal
63 gsl_matrix_char_subrow gsl_matrix_char_subcolumn gsl_matrix_char_view_array
64 gsl_matrix_char_view_array_with_tda gsl_matrix_char_view_vector gsl_matrix_char_view_vector_with_tda
65 gsl_matrix_char_const_submatrix gsl_matrix_char_const_row gsl_matrix_char_const_column
66 gsl_matrix_char_const_diagonal gsl_matrix_char_const_subdiagonal gsl_matrix_char_const_superdiagonal
67 gsl_matrix_char_const_subrow gsl_matrix_char_const_subcolumn gsl_matrix_char_const_view_array
68 gsl_matrix_char_const_view_array_with_tda gsl_matrix_char_const_view_vector gsl_matrix_char_const_view_vector_with_tda
69 gsl_matrix_char_get gsl_matrix_char_set gsl_matrix_char_ptr gsl_matrix_char_const_ptr
70 gsl_matrix_char_set_zero gsl_matrix_char_set_identity
71 gsl_matrix_char_set_all gsl_matrix_char_fread
72 gsl_matrix_char_fwrite gsl_matrix_char_fscanf gsl_matrix_char_fprintf
73 gsl_matrix_char_memcpy gsl_matrix_char_swap
74 gsl_matrix_char_swap_rows gsl_matrix_char_swap_columns
75 gsl_matrix_char_swap_rowcol gsl_matrix_char_transpose gsl_matrix_char_transpose_memcpy
76 gsl_matrix_char_max gsl_matrix_char_min
77 gsl_matrix_char_minmax gsl_matrix_char_max_index
78 gsl_matrix_char_min_index gsl_matrix_char_minmax_index
79 gsl_matrix_char_isnull gsl_matrix_char_ispos gsl_matrix_char_isneg
80 gsl_matrix_char_isnonneg gsl_matrix_char_add
81 gsl_matrix_char_sub gsl_matrix_char_mul_elements gsl_matrix_char_div_elements
82 gsl_matrix_char_scale gsl_matrix_char_add_constant gsl_matrix_char_add_diagonal
83 gsl_matrix_int_alloc gsl_matrix_int_calloc gsl_matrix_int_alloc_from_block
84 gsl_matrix_int_alloc_from_matrix gsl_vector_int_alloc_row_from_matrix gsl_vector_int_alloc_col_from_matrix
85 gsl_matrix_int_free gsl_matrix_int_submatrix gsl_matrix_int_row
86 gsl_matrix_int_column gsl_matrix_int_diagonal gsl_matrix_int_subdiagonal
87 gsl_matrix_int_superdiagonal gsl_matrix_int_subrow gsl_matrix_int_subcolumn gsl_matrix_int_view_array
88 gsl_matrix_int_view_array_with_tda gsl_matrix_int_view_vector gsl_matrix_int_view_vector_with_tda
89 gsl_matrix_int_const_submatrix gsl_matrix_int_const_row gsl_matrix_int_const_column
90 gsl_matrix_int_const_diagonal gsl_matrix_int_const_subdiagonal gsl_matrix_int_const_superdiagonal
91 gsl_matrix_int_const_subrow gsl_matrix_int_const_subcolumn gsl_matrix_int_const_view_array
92 gsl_matrix_int_const_view_array_with_tda gsl_matrix_int_const_view_vector gsl_matrix_int_const_view_vector_with_tda
93 gsl_matrix_int_get gsl_matrix_int_set
94 gsl_matrix_int_ptr gsl_matrix_int_const_ptr
95 gsl_matrix_int_set_zero gsl_matrix_int_set_identity gsl_matrix_int_set_all
96 gsl_matrix_int_fread gsl_matrix_int_fwrite
97 gsl_matrix_int_fscanf gsl_matrix_int_fprintf
98 gsl_matrix_int_memcpy gsl_matrix_int_swap
99 gsl_matrix_int_swap_rows gsl_matrix_int_swap_columns gsl_matrix_int_swap_rowcol
100 gsl_matrix_int_transpose gsl_matrix_int_transpose_memcpy
101 gsl_matrix_int_max gsl_matrix_int_min gsl_matrix_int_minmax
102 gsl_matrix_int_max_index gsl_matrix_int_min_index
103 gsl_matrix_int_minmax_index gsl_matrix_int_isnull
104 gsl_matrix_int_ispos gsl_matrix_int_isneg gsl_matrix_int_isnonneg
105 gsl_matrix_int_add gsl_matrix_int_sub
106 gsl_matrix_int_mul_elements gsl_matrix_int_div_elements gsl_matrix_int_scale
107 gsl_matrix_int_add_constant gsl_matrix_int_add_diagonal
108 gsl_matrix_get_row gsl_matrix_get_col gsl_matrix_set_row gsl_matrix_set_col
109 gsl_matrix_complex_alloc gsl_matrix_complex_calloc gsl_matrix_complex_alloc_from_block
110 gsl_matrix_complex_alloc_from_matrix gsl_vector_complex_alloc_row_from_matrix gsl_vector_complex_alloc_col_from_matrix
111 gsl_matrix_complex_free gsl_matrix_complex_submatrix gsl_matrix_complex_row
112 gsl_matrix_complex_column gsl_matrix_complex_diagonal gsl_matrix_complex_subdiagonal
113 gsl_matrix_complex_superdiagonal gsl_matrix_complex_subrow gsl_matrix_complex_subcolumn
114 gsl_matrix_complex_view_array gsl_matrix_complex_view_array_with_tda gsl_matrix_complex_view_vector
115 gsl_matrix_complex_view_vector_with_tda gsl_matrix_complex_const_submatrix gsl_matrix_complex_const_row
116 gsl_matrix_complex_const_column gsl_matrix_complex_const_diagonal gsl_matrix_complex_const_subdiagonal
117 gsl_matrix_complex_const_superdiagonal gsl_matrix_complex_const_subrow gsl_matrix_complex_const_subcolumn
118 gsl_matrix_complex_const_view_array gsl_matrix_complex_const_view_array_with_tda gsl_matrix_complex_const_view_vector
119 gsl_matrix_complex_const_view_vector_with_tda gsl_matrix_complex_get gsl_matrix_complex_set
120 gsl_matrix_complex_ptr gsl_matrix_complex_const_ptr gsl_matrix_complex_set_zero
121 gsl_matrix_complex_set_identity gsl_matrix_complex_set_all gsl_matrix_complex_fread
122 gsl_matrix_complex_fwrite gsl_matrix_complex_fscanf gsl_matrix_complex_fprintf
123 gsl_matrix_complex_memcpy gsl_matrix_complex_swap gsl_matrix_complex_swap_rows
124 gsl_matrix_complex_swap_columns gsl_matrix_complex_swap_rowcol gsl_matrix_complex_transpose
125 gsl_matrix_complex_transpose_memcpy gsl_matrix_complex_isnull gsl_matrix_complex_ispos
126 gsl_matrix_complex_isneg gsl_matrix_complex_add gsl_matrix_complex_sub
127 gsl_matrix_complex_mul_elements gsl_matrix_complex_div_elements gsl_matrix_complex_scale
128 gsl_matrix_complex_add_constant gsl_matrix_complex_add_diagonal gsl_matrix_complex_get_row
129 gsl_matrix_complex_get_col gsl_matrix_complex_set_row gsl_matrix_complex_set_col /;
133 %EXPORT_TAGS = ( all => [ @EXPORT_OK ],
134 char => [ qw/
135 gsl_matrix_char_alloc
136 gsl_matrix_char_calloc
137 gsl_matrix_char_alloc_from_block
138 gsl_matrix_char_alloc_from_matrix
139 gsl_vector_char_alloc_row_from_matrix
140 gsl_vector_char_alloc_col_from_matrix
141 gsl_matrix_char_free
142 gsl_matrix_char_submatrix
143 gsl_matrix_char_row
144 gsl_matrix_char_column
145 gsl_matrix_char_diagonal
146 gsl_matrix_char_subdiagonal
147 gsl_matrix_char_superdiagonal
148 gsl_matrix_char_subrow
149 gsl_matrix_char_subcolumn
150 gsl_matrix_char_view_array
151 gsl_matrix_char_view_array_with_tda
152 gsl_matrix_char_view_vector
153 gsl_matrix_char_view_vector_with_tda
154 gsl_matrix_char_const_submatrix
155 gsl_matrix_char_const_row
156 gsl_matrix_char_const_column
157 gsl_matrix_char_const_diagonal
158 gsl_matrix_char_const_subdiagonal
159 gsl_matrix_char_const_superdiagonal
160 gsl_matrix_char_const_subrow
161 gsl_matrix_char_const_subcolumn
162 gsl_matrix_char_const_view_array
163 gsl_matrix_char_const_view_array_with_tda
164 gsl_matrix_char_const_view_vector
165 gsl_matrix_char_const_view_vector_with_tda
166 gsl_matrix_char_get
167 gsl_matrix_char_set
168 gsl_matrix_char_ptr
169 gsl_matrix_char_const_ptr
170 gsl_matrix_char_set_zero
171 gsl_matrix_char_set_identity
172 gsl_matrix_char_set_all
173 gsl_matrix_char_fread
174 gsl_matrix_char_fwrite
175 gsl_matrix_char_fscanf
176 gsl_matrix_char_fprintf
177 gsl_matrix_char_memcpy
178 gsl_matrix_char_swap
179 gsl_matrix_char_swap_rows
180 gsl_matrix_char_swap_columns
181 gsl_matrix_char_swap_rowcol
182 gsl_matrix_char_transpose
183 gsl_matrix_char_transpose_memcpy
184 gsl_matrix_char_max
185 gsl_matrix_char_min
186 gsl_matrix_char_minmax
187 gsl_matrix_char_max_index
188 gsl_matrix_char_min_index
189 gsl_matrix_char_minmax_index
190 gsl_matrix_char_isnull
191 gsl_matrix_char_ispos
192 gsl_matrix_char_isneg
193 gsl_matrix_char_isnonneg
194 gsl_matrix_char_add
195 gsl_matrix_char_sub
196 gsl_matrix_char_mul_elements
197 gsl_matrix_char_div_elements
198 gsl_matrix_char_scale
199 gsl_matrix_char_add_constant
200 gsl_matrix_char_add_diagonal
203 double => [ qw/
204 gsl_matrix_alloc
205 gsl_matrix_calloc
206 gsl_matrix_alloc_from_block
207 gsl_matrix_alloc_from_matrix
208 gsl_vector_alloc_row_from_matrix
209 gsl_vector_alloc_col_from_matrix
210 gsl_matrix_free
211 gsl_matrix_submatrix
212 gsl_matrix_row
213 gsl_matrix_column
214 gsl_matrix_diagonal
215 gsl_matrix_subdiagonal
216 gsl_matrix_superdiagonal
217 gsl_matrix_subrow
218 gsl_matrix_subcolumn
219 gsl_matrix_view_array
220 gsl_matrix_view_array_with_tda
221 gsl_matrix_view_vector
222 gsl_matrix_view_vector_with_tda
223 gsl_matrix_const_submatrix
224 gsl_matrix_const_row
225 gsl_matrix_const_column
226 gsl_matrix_const_diagonal
227 gsl_matrix_const_subdiagonal
228 gsl_matrix_const_superdiagonal
229 gsl_matrix_const_subrow
230 gsl_matrix_const_subcolumn
231 gsl_matrix_const_view_array
232 gsl_matrix_const_view_array_with_tda
233 gsl_matrix_const_view_vector
234 gsl_matrix_const_view_vector_with_tda
235 gsl_matrix_get
236 gsl_matrix_set
237 gsl_matrix_ptr
238 gsl_matrix_const_ptr
239 gsl_matrix_set_zero
240 gsl_matrix_set_identity
241 gsl_matrix_set_all
242 gsl_matrix_fread
243 gsl_matrix_fwrite
244 gsl_matrix_fscanf
245 gsl_matrix_fprintf
246 gsl_matrix_memcpy
247 gsl_matrix_swap
248 gsl_matrix_swap_rows
249 gsl_matrix_swap_columns
250 gsl_matrix_swap_rowcol
251 gsl_matrix_transpose
252 gsl_matrix_transpose_memcpy
253 gsl_matrix_max
254 gsl_matrix_minmax
255 gsl_matrix_max_index
256 gsl_matrix_min_index
257 gsl_matrix_minmax_index
258 gsl_matrix_isnull
259 gsl_matrix_ispos
260 gsl_matrix_isneg
261 gsl_matrix_isnonneg
262 gsl_matrix_add
263 gsl_matrix_mul_elements
264 gsl_matrix_div_elements
265 gsl_matrix_scale
266 gsl_matrix_add_constant
267 gsl_matrix_add_diagonal
269 int => [ qw/
270 gsl_matrix_int_alloc
271 gsl_matrix_int_alloc_from_matrix
272 gsl_matrix_int_free
273 gsl_matrix_int_column
274 gsl_matrix_int_superdiagonal
275 gsl_matrix_int_view_array_with_tda
276 gsl_matrix_int_const_submatrix
277 gsl_matrix_int_const_diagonal
278 gsl_matrix_int_const_subrow
279 gsl_matrix_int_const_view_array_with_tda
280 gsl_matrix_int_get
281 gsl_matrix_int_ptr
282 gsl_matrix_int_set_zero
283 gsl_matrix_int_fread
284 gsl_matrix_int_fscanf
285 gsl_matrix_int_memcpy
286 gsl_matrix_int_swap_rows
287 gsl_matrix_int_transpose
288 gsl_matrix_int_max
289 gsl_matrix_int_max_index
290 gsl_matrix_int_minmax_index
291 gsl_matrix_int_ispos
292 gsl_matrix_int_add
293 gsl_matrix_int_mul_elements
294 gsl_matrix_int_add_constant
296 complex => [ qw/
297 gsl_matrix_complex_alloc
298 gsl_matrix_complex_calloc
299 gsl_matrix_complex_alloc_from_block
300 gsl_matrix_complex_alloc_from_matrix
301 gsl_vector_complex_alloc_row_from_matrix
302 gsl_vector_complex_alloc_col_from_matrix
303 gsl_matrix_complex_free
304 gsl_matrix_complex_submatrix
305 gsl_matrix_complex_row
306 gsl_matrix_complex_column
307 gsl_matrix_complex_diagonal
308 gsl_matrix_complex_subdiagonal
309 gsl_matrix_complex_superdiagonal
310 gsl_matrix_complex_subrow
311 gsl_matrix_complex_subcolumn
312 gsl_matrix_complex_view_array
313 gsl_matrix_complex_view_array_with_tda
314 gsl_matrix_complex_view_vector
315 gsl_matrix_complex_view_vector_with_tda
316 gsl_matrix_complex_const_submatrix
317 gsl_matrix_complex_const_row
318 gsl_matrix_complex_const_column
319 gsl_matrix_complex_const_diagonal
320 gsl_matrix_complex_const_subdiagonal
321 gsl_matrix_complex_const_superdiagonal
322 gsl_matrix_complex_const_subrow
323 gsl_matrix_complex_const_subcolumn
324 gsl_matrix_complex_const_view_array
325 gsl_matrix_complex_const_view_array_with_tda
326 gsl_matrix_complex_const_view_vector
327 gsl_matrix_complex_const_view_vector_with_tda
328 gsl_matrix_complex_get
329 gsl_matrix_complex_set
330 gsl_matrix_complex_ptr
331 gsl_matrix_complex_const_ptr
332 gsl_matrix_complex_set_zero gsl_matrix_complex_set_identity
333 gsl_matrix_complex_set_all
334 gsl_matrix_complex_fread
335 gsl_matrix_complex_fwrite
336 gsl_matrix_complex_fscanf
337 gsl_matrix_complex_fprintf
338 gsl_matrix_complex_memcpygsl_matrix_complex_swap
339 gsl_matrix_complex_swap_rows
340 gsl_matrix_complex_swap_columns
341 gsl_matrix_complex_swap_rowcol
342 gsl_matrix_complex_transpose
343 gsl_matrix_complex_transpose_memcpy
344 gsl_matrix_complex_isnull
345 gsl_matrix_complex_ispos
346 gsl_matrix_complex_isneg
347 gsl_matrix_complex_add
348 gsl_matrix_complex_sub
349 gsl_matrix_complex_mul_elements
350 gsl_matrix_complex_div_elements
351 gsl_matrix_complex_scale
352 gsl_matrix_complex_add_constant
353 gsl_matrix_complex_add_diagonal
354 gsl_matrix_complex_get_row
355 gsl_matrix_complex_get_col
356 gsl_matrix_complex_set_row
357 gsl_matrix_complex_set_col
358 /]);
360 =head1 NAME
362 Math::GSL::Matrix - Mathematical functions concerning Matrices
364 =head1 SYNOPSIS
366 use Math::GSL::Matrix qw/:all/;
367 my $matrix1 = Math::GSL::Matrix->new(5,5); # OO interface
368 my $matrix2 = gsl_matrix_alloc(5,5); # standard interface
371 =head1 Objected Oriented Interface to GSL Math::GSL::Matrix
373 =head2 Math::GSL::Matrix->new()
375 Creates a new Matrix of the given size.
377 my $matrix = Math::GSL::Matrix->new(10,10);
378 =cut
380 sub new
382 my ($class, $rows, $cols) = @_;
383 my $this = {};
384 my $matrix;
385 if ( defined $rows && defined $cols &&
386 $rows > 0 && $cols > 0 &&
387 (int $rows == $rows) && (int $cols == $cols)){
389 $matrix = gsl_matrix_alloc($rows,$cols);
390 } else {
391 croak( __PACKAGE__.'::new($x,$y) - $x and $y must be positive integers');
393 gsl_matrix_set_zero($matrix);
394 $this->{_matrix} = $matrix;
395 ($this->{_rows}, $this->{_cols}) = ($rows,$cols);
396 bless $this, $class;
398 =head2 raw()
400 Get the underlying GSL matrix object created by SWIG, useful for using gsl_matrix_* functions which do not have an OO counterpart.
402 my $matrix = Math::GSL::Matrix->new(3,3);
403 my $gsl_matrix = $matrix->raw;
404 my $stuff = gsl_matrix_get($gsl_matrix, 1, 2);
406 =cut
407 sub raw { (shift)->{_matrix} }
408 =head2 rows()
410 Returns the number of rows in the matrix.
412 my $rows = $matrix->rows;
413 =cut
415 sub rows { (shift)->{_rows} }
417 =head2 cols()
419 Returns the number of columns in the matrix.
421 my $cols = $matrix->cols;
422 =cut
424 sub cols { (shift)->{_cols} }
426 =head2 as_list()
428 Get the contents of a Math::GSL::Matrix object as a Perl list.
430 my $matrix = Math::GSL::Matrix->new(3,3);
432 my @matrix = $matrix->as_list;
433 =cut
436 sub as_list
438 my $self = shift;
439 my $line;
440 my @part;
441 my @total;
442 for($line=0; $line<$self->rows; $line++){
443 @part = map {
444 gsl_matrix_get($self->raw, $line, $_)
445 } (0 .. $self->cols-1 );
446 push(@total, @part);
448 return @total;
451 =head2 row()
453 Returns a row matrix of the row you enter.
455 my $matrix = Math::GSL::Matrix->new(3,3);
457 my $matrix_row = $matrix->row(0);
459 =cut
461 sub row
463 my ($self, $row) = @_;
464 croak (__PACKAGE__.'::$matrix->row($row) - invalid $row value')
465 unless (($row < $self->rows) and $row >= 0);
467 my $rowmat = Math::GSL::Matrix->new(1,$self->cols);
469 my @got;
470 for my $n (0 .. $self->cols-1) {
471 push (@got, gsl_matrix_get($self->raw, $row, $n));
473 for my $n (0 .. $self->cols-1) {
474 gsl_matrix_set($rowmat->raw, 0, $n, $got[$n]); }
476 return $rowmat;
479 =head2 col()
481 Returns a col matrix of the column you enter.
483 my $matrix = Math::GSL::Matrix->new(3,3);
485 my $matrix_col = $matrix->col(0);
487 =cut
489 sub col
491 my ($self, $col) = @_;
492 croak (__PACKAGE__."::\$matrix->col(\$col) - $col not a valid column")
493 unless ($col < $self->cols and $col >= 0);
495 my $colvec = Math::GSL::Vector->new($self->cols);
496 my $colmat = Math::GSL::Matrix->new($self->rows, 1);
498 my @got;
499 for my $n (0 .. $self->rows-1) {
500 push (@got, gsl_matrix_get($self->raw, $n, $col));
502 for my $n (0 .. $self->rows-1) {
503 gsl_matrix_set($colmat->raw, $n, 0, $got[$n]); }
505 return $colmat;
508 =head2 set_row()
510 Sets a the values of a row with the elements of an array.
512 my $matrix = Math::GSL::Matrix->new(3,3);
513 $matrix->set_row(0, [8, 6, 2]);
515 You can also set multiple rows at once with chained calls:
516 my $matrix = Math::GSL::Matrix->new(3,3);
517 $matrix->set_row(0, [8, 6, 2])
518 ->set_row(1, [2, 4, 1]);
521 =cut
523 sub set_row {
524 my ($self, $row, $values) = @_;
525 my $length = $#$values;
526 die __PACKAGE__.'::new($x, $values) - $values must be a nonempty array reference' if $length == -1;
527 die __PACKAGE__.'::set_row($x, $values) - $x must be a valid row number' if ($row < 0 || $row >= $self->rows);
528 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);
529 map { gsl_matrix_set($self->raw, $row, $_, $values->[$_]) } (0..$length);
530 return $self;
533 =head2 set_col()
535 Sets a the values of a column with the elements of an array.
537 my $matrix = Math::GSL::Matrix->new(3,3);
538 $matrix->set_col(0, [8, 6, 2]);
540 You can also set multiple columns at once with chained calls:
541 my $matrix = Math::GSL::Matrix->new(3,3);
542 $matrix->set_col(0, [8, 6, 2])
543 ->set_col(1, [2, 4, 1]);
546 =cut
548 sub set_col {
549 my ($self, $col, $values) = @_;
550 my $length = $#$values;
551 die __PACKAGE__.'::new($x, $values) - $values must be a nonempty array reference' if $length == -1;
552 die __PACKAGE__.'::set_row($x, $values) - $x must be a valid column number' if ($col < 0 || $col >= $self->cols);
553 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);
554 map { gsl_matrix_set($self->raw, $_, $col, $values->[$_]) } (0..$length);
555 return $self;
558 =head1 DESCRIPTION
560 Here is a list of all the functions included in this module :
562 =over 1
564 =item C<gsl_matrix_alloc($i, $j)> - Return a gsl_matrix of $i rows and $j columns
566 =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
568 =item C<gsl_matrix_alloc_from_block> -
570 =item C<gsl_matrix_free> -
572 =item C<gsl_matrix_alloc_from_matrix > -
574 =item C<gsl_vector_alloc_row_from_matrix> -
576 =item C<gsl_vector_alloc_col_from_matrix > -
578 =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.
580 =item C<gsl_matrix_row($m , $i)> - Return a vector view of the $i-th row of the matrix $m
582 =item C<gsl_matrix_column($m, $j)> - Return a vector view of the $j-th column of the matrix $m
584 =item C<gsl_matrix_diagonal($m)> - Return a vector view of the diagonal of the vector. The matrix doesn't have to be square.
586 =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.
588 =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.
590 =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.
592 =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.
594 =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.
596 =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.
598 =item C<gsl_matrix_view_vector> -
600 =item C<gsl_matrix_view_vector_with_tda> -
602 =item C<gsl_matrix_const_submatrix> -
604 =item C<gsl_matrix_get($m, $i, $j)> - Return the (i,j)-th element of the matrix $m
606 =item C<gsl_matrix_set($m, $i, $j, $x)> - Set the value of the (i,j)-th element of the matrix $m to $x
608 =item C<gsl_matrix_ptr> -
610 =item C<gsl_matrix_const_ptr> -
612 =item C<gsl_matrix_set_zero($m)> - Set all the elements of the matrix $m to zero
614 =item C<gsl_matrix_set_identity($m)> - Set the elements of the matrix $m to the corresponding elements of the identity matrix
616 =item C<gsl_matrix_set_all($m, $x)> - Set all the elements of the matrix $m to the value $x
618 =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
620 =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
622 =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
624 =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
626 =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.
628 =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.
630 =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.
632 =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.
634 =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.
636 =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.
638 =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.
640 =item C<gsl_matrix_max($m)> - Return the maximum value in the matrix $m
642 =item C<gsl_matrix_min($m)> - Return the minimum value in the matrix $m
644 =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.
646 =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.
648 =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.
650 =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
652 =item C<gsl_matrix_isnull($m)> - Return 1 if all the elements of the matrix $m are zero, 0 otherwise
654 =item C<gsl_matrix_ispos($m)> - Return 1 if all the elements of the matrix $m are strictly positve, 0 otherwise
656 =item C<gsl_matrix_isneg($m)> - Return 1 if all the elements of the matrix $m are strictly negative, 0 otherwise
658 =item C<gsl_matrix_isnonneg($m)> - Return 1 if all the elements of the matrix $m are non-negatuive, 0 otherwise
660 =item C<gsl_matrix_add($a, $b)> - Add the elements of matrix $b to the elements of matrix $a
662 =item C<gsl_matrix_sub($a, $b)> - Subtract the elements of matrix $b from the elements of matrix $a
664 =item C<gsl_matrix_mul_elements($a, $b)> - Multiplie the elements of matrix $a by the elements of matrix $b
666 =item C<gsl_matrix_div_elements($a, $b)> - Divide the elements of matrix $a by the elements of matrix $b
668 =item C<gsl_matrix_scale($a, $x)> - Multiplie the elements of matrix $a by the constant factor $x
670 =item C<gsl_matrix_add_constant($a, $x)> - Add the constant value $x to the elements of the matrix $a
672 =item C<gsl_matrix_add_diagonal($a, $x)> - Add the constant value $x to the elements of the diagonal of the matrix $a
674 =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.
676 =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.
678 =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.
680 =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.
682 =back
684 The following functions are specific to matrices containing complex numbers :
686 =over 1
688 =item C<gsl_matrix_complex_alloc >
690 =item C<gsl_matrix_complex_calloc >
692 =item C<gsl_matrix_complex_alloc_from_block >
694 =item C<gsl_matrix_complex_alloc_from_matrix >
696 =item C<gsl_vector_complex_alloc_row_from_matrix >
698 =item C<gsl_vector_complex_alloc_col_from_matrix >
700 =item C<gsl_matrix_complex_free >
702 =item C<gsl_matrix_complex_submatrix >
704 =item C<gsl_matrix_complex_row >
706 =item C<gsl_matrix_complex_column >
708 =item C<gsl_matrix_complex_diagonal >
710 =item C<gsl_matrix_complex_subdiagonal >
712 =item C<gsl_matrix_complex_superdiagonal >
714 =item C<gsl_matrix_complex_subrow >
716 =item C<gsl_matrix_complex_subcolumn >
718 =item C<gsl_matrix_complex_view_array >
720 =item C<gsl_matrix_complex_view_array_with_tda >
722 =item C<gsl_matrix_complex_view_vector >
724 =item C<gsl_matrix_complex_view_vector_with_tda >
726 =item C<gsl_matrix_complex_const_submatrix >
728 =item C<gsl_matrix_complex_const_row >
730 =item C<gsl_matrix_complex_const_column >
732 =item C<gsl_matrix_complex_const_diagonal >
734 =item C<gsl_matrix_complex_const_subdiagonal >
736 =item C<gsl_matrix_complex_const_superdiagonal >
738 =item C<gsl_matrix_complex_const_subrow >
740 =item C<gsl_matrix_complex_const_subcolumn >
742 =item C<gsl_matrix_complex_const_view_array >
744 =item C<gsl_matrix_complex_const_view_array_with_tda >
746 =item C<gsl_matrix_complex_const_view_vector >
748 =item C<gsl_matrix_complex_const_view_vector_with_tda >
750 =item C<gsl_matrix_complex_get>
752 =item C<gsl_matrix_complex_set>
754 =item C<gsl_matrix_complex_ptr>
756 =item C<gsl_matrix_complex_const_ptr>
758 =item C<gsl_matrix_complex_set_zero >
760 =item C<gsl_matrix_complex_set_identity >
762 =item C<gsl_matrix_complex_set_all >
764 =item C<gsl_matrix_complex_fread >
766 =item C<gsl_matrix_complex_fwrite >
768 =item C<gsl_matrix_complex_fscanf >
770 =item C<gsl_matrix_complex_fprintf >
772 =item C<gsl_matrix_complex_memcpy>
774 =item C<gsl_matrix_complex_swap>
776 =item C<gsl_matrix_complex_swap_rows>
778 =item C<gsl_matrix_complex_swap_columns>
780 =item C<gsl_matrix_complex_swap_rowcol>
782 =item C<gsl_matrix_complex_transpose >
784 =item C<gsl_matrix_complex_transpose_memcpy >
786 =item C<gsl_matrix_complex_isnull >
788 =item C<gsl_matrix_complex_ispos >
790 =item C<gsl_matrix_complex_isneg >
792 =item C<gsl_matrix_complex_add >
794 =item C<gsl_matrix_complex_sub >
796 =item C<gsl_matrix_complex_mul_elements >
798 =item C<gsl_matrix_complex_div_elements >
800 =item C<gsl_matrix_complex_scale >
802 =item C<gsl_matrix_complex_add_constant >
804 =item C<gsl_matrix_complex_add_diagonal >
806 =item C<gsl_matrix_complex_get_row>
808 =item C<gsl_matrix_complex_get_col>
810 =item C<gsl_matrix_complex_set_row>
812 =item C<gsl_matrix_complex_set_col>
814 =back
816 These are related to constant views of a matrix.
818 =over 1
821 =item C<gsl_matrix_const_row>
823 =item C<gsl_matrix_const_colum>
825 =item C<gsl_matrix_const_diagonal>
827 =item C<gsl_matrix_const_subdiagonal>
829 =item C<gsl_matrix_const_superdiagonal>
831 =item C<gsl_matrix_const_subrow>
833 =item C<gsl_matrix_const_subcolumn>
835 =item C<gsl_matrix_const_view_array>
837 =item C<gsl_matrix_const_view_array_with_tda>
839 =back
842 The following functions are similar to those above but work with C<char>'s and C<int>'s. We are not quite
843 sure if anyone wants these. Please speak up if you do and/or submit some patches to this documentation, please!
845 =over 1
848 =item gsl_matrix_const_view_vector
850 =item gsl_matrix_const_view_vector_with_tda
852 =item gsl_matrix_char_alloc
854 =item gsl_matrix_char_calloc
856 =item gsl_matrix_char_alloc_from_block
858 =item gsl_matrix_char_alloc_from_matrix
860 =item gsl_vector_char_alloc_row_from_matrix
862 =item gsl_vector_char_alloc_col_from_matrix
864 =item gsl_matrix_char_free
866 =item gsl_matrix_char_submatrix
868 =item gsl_matrix_char_row
870 =item gsl_matrix_char_column
872 =item gsl_matrix_char_diagonal
874 =item gsl_matrix_char_subdiagonal
876 =item gsl_matrix_char_superdiagonal
878 =item gsl_matrix_char_subrow
880 =item gsl_matrix_char_subcolumn
882 =item gsl_matrix_char_view_array
884 =item gsl_matrix_char_view_array_with_tda
886 =item gsl_matrix_char_view_vector
888 =item gsl_matrix_char_view_vector_with_tda
890 =item gsl_matrix_char_const_submatrix
892 =item gsl_matrix_char_const_row
894 =item gsl_matrix_char_const_column
896 =item gsl_matrix_char_const_diagonal
898 =item gsl_matrix_char_const_subdiagonal
900 =item gsl_matrix_char_const_superdiagonal
902 =item gsl_matrix_char_const_subrow
904 =item gsl_matrix_char_const_subcolumn
906 =item gsl_matrix_char_const_view_array
908 =item gsl_matrix_char_const_view_array_with_tda
910 =item gsl_matrix_char_const_view_vector
912 =item gsl_matrix_char_const_view_vector_with_tda
914 =item gsl_matrix_char_get
916 =item gsl_matrix_char_set
918 =item gsl_matrix_char_ptr
920 =item gsl_matrix_char_const_ptr
922 =item gsl_matrix_char_set_zero
924 =item gsl_matrix_char_set_identity
926 =item gsl_matrix_char_set_all
928 =item gsl_matrix_char_fread
930 =item gsl_matrix_char_fwrite
932 =item gsl_matrix_char_fscanf
934 =item gsl_matrix_char_fprintf
936 =item gsl_matrix_char_memcpy
938 =item gsl_matrix_char_swap
940 =item gsl_matrix_char_swap_rows
942 =item gsl_matrix_char_swap_columns
944 =item gsl_matrix_char_swap_rowcol
946 =item gsl_matrix_char_transpose
948 =item gsl_matrix_char_transpose_memcpy
950 =item gsl_matrix_char_max
952 =item gsl_matrix_char_min
954 =item gsl_matrix_char_minmax
956 =item gsl_matrix_char_max_index
958 =item gsl_matrix_char_min_index
960 =item gsl_matrix_char_minmax_index
962 =item gsl_matrix_char_isnull
964 =item gsl_matrix_char_ispos
966 =item gsl_matrix_char_isneg
968 =item gsl_matrix_char_isnonneg
970 =item gsl_matrix_char_add
972 =item gsl_matrix_char_sub
974 =item gsl_matrix_char_mul_elements
976 =item gsl_matrix_char_div_elements
978 =item gsl_matrix_char_scale
980 =item gsl_matrix_char_add_constant
982 =item gsl_matrix_char_add_diagonal
984 =item gsl_matrix_int_alloc
986 =item gsl_matrix_int_calloc
988 =item gsl_matrix_int_alloc_from_block
990 =item gsl_matrix_int_alloc_from_matrix
992 =item gsl_vector_int_alloc_row_from_matrix
994 =item gsl_vector_int_alloc_col_from_matrix
996 =item gsl_matrix_int_free
998 =item gsl_matrix_int_submatrix
1000 =item gsl_matrix_int_row
1002 =item gsl_matrix_int_column
1004 =item gsl_matrix_int_diagonal
1006 =item gsl_matrix_int_subdiagonal
1008 =item gsl_matrix_int_superdiagonal
1010 =item gsl_matrix_int_subrow
1012 =item gsl_matrix_int_subcolumn
1014 =item gsl_matrix_int_view_array
1016 =item gsl_matrix_int_view_array_with_tda
1018 =item gsl_matrix_int_view_vector
1020 =item gsl_matrix_int_view_vector_with_tda
1022 =item gsl_matrix_int_const_submatrix
1024 =item gsl_matrix_int_const_row
1026 =item gsl_matrix_int_const_column
1028 =item gsl_matrix_int_ptr
1030 =item gsl_matrix_int_const_ptr
1032 =item gsl_matrix_int_set_zero
1034 =item gsl_matrix_int_set_identity
1036 =item gsl_matrix_int_set_all
1038 =item gsl_matrix_int_fread
1040 =item gsl_matrix_int_fwrite
1042 =item gsl_matrix_int_fscanf
1044 =item gsl_matrix_int_fprintf
1046 =item gsl_matrix_int_memcpy
1048 =item gsl_matrix_int_swap
1050 =item gsl_matrix_int_swap_rows
1052 =item gsl_matrix_int_swap_columns
1054 =item gsl_matrix_int_swap_rowcol
1056 =item gsl_matrix_int_transpose
1058 =item gsl_matrix_int_transpose_memcpy
1060 =item gsl_matrix_int_max
1062 =item gsl_matrix_int_min
1064 =item gsl_matrix_int_minmax
1066 =item gsl_matrix_int_max_index
1068 =item gsl_matrix_int_min_index
1070 =item gsl_matrix_int_minmax_index
1072 =item gsl_matrix_int_isnull
1074 =item gsl_matrix_int_ispos
1076 =item gsl_matrix_int_isneg
1078 =item gsl_matrix_int_isnonneg
1080 =item gsl_matrix_int_add
1082 =item gsl_matrix_int_sub
1084 =item gsl_matrix_int_mul_elements
1086 =item gsl_matrix_int_div_elements
1088 =item gsl_matrix_int_scale
1090 =item gsl_matrix_int_add_constant
1092 =item gsl_matrix_int_add_diagonal
1094 =back
1096 You have to add the functions you want to use inside the qw /put_funtion_here /.
1097 You can also write use Math::GSL::Matrix qw/:all/ to use all avaible functions of the module.
1098 Other tags are also avaible, here is a complete list of all tags for this module :
1100 =over 1
1102 =item C<all>
1104 =item C<int>
1106 =item C<double>
1108 =item C<char>
1110 =item C<complex>
1112 =back
1114 For more informations on the functions, we refer you to the GSL offcial documentation
1115 L<http://www.gnu.org/software/gsl/manual/html_node/>
1117 Tip : search on google: site:http://www.gnu.org/software/gsl/manual/html_node/name_of_the_function_you_want
1119 =head1 EXAMPLES
1121 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>
1123 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.
1125 use Math::GSL::Matrix qw/:all/;
1126 my $m = gsl_matrix_alloc (10,3);
1127 for my $i (0..9){
1128 for my $j (0..2){
1129 gsl_matrix_set($m, $i, $j, 0.23 + 100*$i + $j);
1133 for my $i (0..99){ # OUT OF RANGE ERROR
1134 for my $j (0..2){
1135 print "m($i, $j) = " . gsl_matrix_get ($m, $i, $j) . "\n";
1138 gsl_matrix_free ($m);
1141 use Math::GSL::Matrix qw/:all/;
1143 my $m = gsl_matrix_alloc (100, 100);
1144 my $a = gsl_matrix_alloc (100, 100);
1146 for my $i (0..99){
1147 for my $j (0..99){
1148 gsl_matrix_set ($m, $i, $j, 0.23 + $i + $j);
1152 The next program shows how to write a matrix to a file.
1154 my $out = gsl_fopen("test.dat", "wb");
1155 gsl_matrix_fwrite ($out, $m);
1156 gsl_fclose ($out);
1158 my $in = gsl_fopen("test.dat", "rb");
1159 gsl_matrix_fread ($in, $a);
1160 gsl_fclose($in);
1162 my $k=0;
1163 for my $i (0..99){
1164 for my $j (0..99){
1165 $mij = gsl_matrix_get ($m, $i, $j);
1166 $aij = gsl_matrix_get ($a, $i, $j);
1167 $k++ if ($mij != $aij);
1171 gsl_matrix_free($m);
1172 gsl_matrix_free($a);
1174 print "differences = $k (should be zero)\n";
1176 =head1 AUTHORS
1178 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
1180 =head1 COPYRIGHT AND LICENSE
1182 Copyright (C) 2008 Jonathan Leto and Thierry Moisan
1184 This program is free software; you can redistribute it and/or modify it
1185 under the same terms as Perl itself.
1187 =cut