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