Adding documentation for the view_array function of Vector
[Math-GSL.git] / Matrix.i
blob83ca3758413bb7f852447e2f04b8884de08afe97
1 %module "Math::GSL::Matrix"
3 %include "typemaps.i"
4 %include "gsl_typemaps.i"
6 %apply int *OUTPUT { size_t *imin, size_t *imax, size_t *jmin, size_t *jmax };
7 %apply double *OUTPUT { double * min_out, double * max_out };
9 FILE * fopen(char *, char *);
10 int fclose(FILE *);
13 #include "gsl/gsl_matrix.h"
14 #include "gsl/gsl_complex.h"
15 #include "gsl/gsl_vector_double.h"
16 #include "gsl/gsl_matrix_double.h"
17 #include "gsl/gsl_matrix_int.h"
18 #include "gsl/gsl_matrix_complex_double.h"
19 #include "gsl/gsl_matrix_char.h"
22 %include "gsl/gsl_matrix.h"
23 %include "gsl/gsl_complex.h"
24 %include "gsl/gsl_vector_double.h"
25 %include "gsl/gsl_matrix_double.h"
26 %include "gsl/gsl_matrix_int.h"
27 %include "gsl/gsl_matrix_complex_double.h"
28 %include "gsl/gsl_matrix_char.h"
31 %perlcode %{
33 no warnings 'redefine';
34 use Carp qw/croak/;
35 use Math::GSL qw/:all/;
36 use Math::GSL::Errno qw/:all/;
38 @EXPORT_OK = qw/fopen fclose
39 gsl_matrix_alloc gsl_matrix_calloc gsl_matrix_alloc_from_block
40 gsl_matrix_alloc_from_matrix gsl_vector_alloc_row_from_matrix
41 gsl_vector_alloc_col_from_matrix gsl_matrix_free gsl_matrix_submatrix
42 gsl_matrix_row gsl_matrix_column gsl_matrix_diagonal
43 gsl_matrix_subdiagonal gsl_matrix_superdiagonal gsl_matrix_subrow
44 gsl_matrix_subcolumn gsl_matrix_view_array
45 gsl_matrix_view_array_with_tda gsl_matrix_view_vector
46 gsl_matrix_view_vector_with_tda gsl_matrix_const_submatrix
47 gsl_matrix_const_row gsl_matrix_const_column gsl_matrix_const_diagonal
48 gsl_matrix_const_subdiagonal gsl_matrix_const_superdiagonal
49 gsl_matrix_const_subrow gsl_matrix_const_subcolumn
50 gsl_matrix_const_view_array gsl_matrix_const_view_array_with_tda
51 gsl_matrix_const_view_vector gsl_matrix_const_view_vector_with_tda
52 gsl_matrix_get gsl_matrix_set gsl_matrix_ptr gsl_matrix_const_ptr
53 gsl_matrix_set_zero gsl_matrix_set_identity gsl_matrix_set_all
54 gsl_matrix_fread gsl_matrix_fwrite gsl_matrix_fscanf gsl_matrix_fprintf
55 gsl_matrix_memcpy gsl_matrix_swap gsl_matrix_swap_rows
56 gsl_matrix_swap_columns gsl_matrix_swap_rowcol gsl_matrix_transpose
57 gsl_matrix_transpose_memcpy gsl_matrix_max gsl_matrix_min
58 gsl_matrix_minmax gsl_matrix_max_index gsl_matrix_min_index
59 gsl_matrix_minmax_index gsl_matrix_isnull gsl_matrix_ispos
60 gsl_matrix_isneg gsl_matrix_isnonneg gsl_matrix_add gsl_matrix_sub
61 gsl_matrix_mul_elements gsl_matrix_div_elements gsl_matrix_scale
62 gsl_matrix_add_constant gsl_matrix_add_diagonal
63 gsl_matrix_char_alloc gsl_matrix_char_calloc gsl_matrix_char_alloc_from_block
64 gsl_matrix_char_alloc_from_matrix gsl_vector_char_alloc_row_from_matrix gsl_vector_char_alloc_col_from_matrix
65 gsl_matrix_char_free gsl_matrix_char_submatrix
66 gsl_matrix_char_row gsl_matrix_char_column
67 gsl_matrix_char_diagonal gsl_matrix_char_subdiagonal gsl_matrix_char_superdiagonal
68 gsl_matrix_char_subrow gsl_matrix_char_subcolumn gsl_matrix_char_view_array
69 gsl_matrix_char_view_array_with_tda gsl_matrix_char_view_vector gsl_matrix_char_view_vector_with_tda
70 gsl_matrix_char_const_submatrix gsl_matrix_char_const_row gsl_matrix_char_const_column
71 gsl_matrix_char_const_diagonal gsl_matrix_char_const_subdiagonal gsl_matrix_char_const_superdiagonal
72 gsl_matrix_char_const_subrow gsl_matrix_char_const_subcolumn gsl_matrix_char_const_view_array
73 gsl_matrix_char_const_view_array_with_tda gsl_matrix_char_const_view_vector gsl_matrix_char_const_view_vector_with_tda
74 gsl_matrix_char_get gsl_matrix_char_set gsl_matrix_char_ptr gsl_matrix_char_const_ptr
75 gsl_matrix_char_set_zero gsl_matrix_char_set_identity
76 gsl_matrix_char_set_all gsl_matrix_char_fread
77 gsl_matrix_char_fwrite gsl_matrix_char_fscanf gsl_matrix_char_fprintf
78 gsl_matrix_char_memcpy gsl_matrix_char_swap
79 gsl_matrix_char_swap_rows gsl_matrix_char_swap_columns
80 gsl_matrix_char_swap_rowcol gsl_matrix_char_transpose gsl_matrix_char_transpose_memcpy
81 gsl_matrix_char_max gsl_matrix_char_min
82 gsl_matrix_char_minmax gsl_matrix_char_max_index
83 gsl_matrix_char_min_index gsl_matrix_char_minmax_index
84 gsl_matrix_char_isnull gsl_matrix_char_ispos gsl_matrix_char_isneg
85 gsl_matrix_char_isnonneg gsl_matrix_char_add
86 gsl_matrix_char_sub gsl_matrix_char_mul_elements gsl_matrix_char_div_elements
87 gsl_matrix_char_scale gsl_matrix_char_add_constant gsl_matrix_char_add_diagonal
88 gsl_matrix_int_alloc gsl_matrix_int_calloc gsl_matrix_int_alloc_from_block
89 gsl_matrix_int_alloc_from_matrix gsl_vector_int_alloc_row_from_matrix gsl_vector_int_alloc_col_from_matrix
90 gsl_matrix_int_free gsl_matrix_int_submatrix gsl_matrix_int_row
91 gsl_matrix_int_column gsl_matrix_int_diagonal gsl_matrix_int_subdiagonal
92 gsl_matrix_int_superdiagonal gsl_matrix_int_subrow gsl_matrix_int_subcolumn gsl_matrix_int_view_array
93 gsl_matrix_int_view_array_with_tda gsl_matrix_int_view_vector gsl_matrix_int_view_vector_with_tda
94 gsl_matrix_int_const_submatrix gsl_matrix_int_const_row gsl_matrix_int_const_column
95 gsl_matrix_int_const_diagonal gsl_matrix_int_const_subdiagonal gsl_matrix_int_const_superdiagonal
96 gsl_matrix_int_const_subrow gsl_matrix_int_const_subcolumn gsl_matrix_int_const_view_array
97 gsl_matrix_int_const_view_array_with_tda gsl_matrix_int_const_view_vector gsl_matrix_int_const_view_vector_with_tda
98 gsl_matrix_int_get gsl_matrix_int_set
99 gsl_matrix_int_ptr gsl_matrix_int_const_ptr
100 gsl_matrix_int_set_zero gsl_matrix_int_set_identity gsl_matrix_int_set_all
101 gsl_matrix_int_fread gsl_matrix_int_fwrite
102 gsl_matrix_int_fscanf gsl_matrix_int_fprintf
103 gsl_matrix_int_memcpy gsl_matrix_int_swap
104 gsl_matrix_int_swap_rows gsl_matrix_int_swap_columns gsl_matrix_int_swap_rowcol
105 gsl_matrix_int_transpose gsl_matrix_int_transpose_memcpy
106 gsl_matrix_int_max gsl_matrix_int_min gsl_matrix_int_minmax
107 gsl_matrix_int_max_index gsl_matrix_int_min_index
108 gsl_matrix_int_minmax_index gsl_matrix_int_isnull
109 gsl_matrix_int_ispos gsl_matrix_int_isneg gsl_matrix_int_isnonneg
110 gsl_matrix_int_add gsl_matrix_int_sub
111 gsl_matrix_int_mul_elements gsl_matrix_int_div_elements gsl_matrix_int_scale
112 gsl_matrix_int_add_constant gsl_matrix_int_add_diagonal
113 gsl_matrix_get_row gsl_matrix_get_col gsl_matrix_set_row gsl_matrix_set_col
114 gsl_matrix_complex_alloc gsl_matrix_complex_calloc gsl_matrix_complex_alloc_from_block
115 gsl_matrix_complex_alloc_from_matrix gsl_vector_complex_alloc_row_from_matrix gsl_vector_complex_alloc_col_from_matrix
116 gsl_matrix_complex_free gsl_matrix_complex_submatrix gsl_matrix_complex_row
117 gsl_matrix_complex_column gsl_matrix_complex_diagonal gsl_matrix_complex_subdiagonal
118 gsl_matrix_complex_superdiagonal gsl_matrix_complex_subrow gsl_matrix_complex_subcolumn
119 gsl_matrix_complex_view_array gsl_matrix_complex_view_array_with_tda gsl_matrix_complex_view_vector
120 gsl_matrix_complex_view_vector_with_tda gsl_matrix_complex_const_submatrix gsl_matrix_complex_const_row
121 gsl_matrix_complex_const_column gsl_matrix_complex_const_diagonal gsl_matrix_complex_const_subdiagonal
122 gsl_matrix_complex_const_superdiagonal gsl_matrix_complex_const_subrow gsl_matrix_complex_const_subcolumn
123 gsl_matrix_complex_const_view_array gsl_matrix_complex_const_view_array_with_tda gsl_matrix_complex_const_view_vector
124 gsl_matrix_complex_const_view_vector_with_tda gsl_matrix_complex_get gsl_matrix_complex_set
125 gsl_matrix_complex_ptr gsl_matrix_complex_const_ptr gsl_matrix_complex_set_zero
126 gsl_matrix_complex_set_identity gsl_matrix_complex_set_all gsl_matrix_complex_fread
127 gsl_matrix_complex_fwrite gsl_matrix_complex_fscanf gsl_matrix_complex_fprintf
128 gsl_matrix_complex_memcpy gsl_matrix_complex_swap gsl_matrix_complex_swap_rows
129 gsl_matrix_complex_swap_columns gsl_matrix_complex_swap_rowcol gsl_matrix_complex_transpose
130 gsl_matrix_complex_transpose_memcpy gsl_matrix_complex_isnull gsl_matrix_complex_ispos
131 gsl_matrix_complex_isneg gsl_matrix_complex_add gsl_matrix_complex_sub
132 gsl_matrix_complex_mul_elements gsl_matrix_complex_div_elements gsl_matrix_complex_scale
133 gsl_matrix_complex_add_constant gsl_matrix_complex_add_diagonal gsl_matrix_complex_get_row
134 gsl_matrix_complex_get_col gsl_matrix_complex_set_row gsl_matrix_complex_set_col /;
138 %EXPORT_TAGS = ( all => [ @EXPORT_OK ],
139 char => [ qw/
140 gsl_matrix_char_alloc
141 gsl_matrix_char_calloc
142 gsl_matrix_char_alloc_from_block
143 gsl_matrix_char_alloc_from_matrix
144 gsl_vector_char_alloc_row_from_matrix
145 gsl_vector_char_alloc_col_from_matrix
146 gsl_matrix_char_free
147 gsl_matrix_char_submatrix
148 gsl_matrix_char_row
149 gsl_matrix_char_column
150 gsl_matrix_char_diagonal
151 gsl_matrix_char_subdiagonal
152 gsl_matrix_char_superdiagonal
153 gsl_matrix_char_subrow
154 gsl_matrix_char_subcolumn
155 gsl_matrix_char_view_array
156 gsl_matrix_char_view_array_with_tda
157 gsl_matrix_char_view_vector
158 gsl_matrix_char_view_vector_with_tda
159 gsl_matrix_char_const_submatrix
160 gsl_matrix_char_const_row
161 gsl_matrix_char_const_column
162 gsl_matrix_char_const_diagonal
163 gsl_matrix_char_const_subdiagonal
164 gsl_matrix_char_const_superdiagonal
165 gsl_matrix_char_const_subrow
166 gsl_matrix_char_const_subcolumn
167 gsl_matrix_char_const_view_array
168 gsl_matrix_char_const_view_array_with_tda
169 gsl_matrix_char_const_view_vector
170 gsl_matrix_char_const_view_vector_with_tda
171 gsl_matrix_char_get
172 gsl_matrix_char_set
173 gsl_matrix_char_ptr
174 gsl_matrix_char_const_ptr
175 gsl_matrix_char_set_zero
176 gsl_matrix_char_set_identity
177 gsl_matrix_char_set_all
178 gsl_matrix_char_fread
179 gsl_matrix_char_fwrite
180 gsl_matrix_char_fscanf
181 gsl_matrix_char_fprintf
182 gsl_matrix_char_memcpy
183 gsl_matrix_char_swap
184 gsl_matrix_char_swap_rows
185 gsl_matrix_char_swap_columns
186 gsl_matrix_char_swap_rowcol
187 gsl_matrix_char_transpose
188 gsl_matrix_char_transpose_memcpy
189 gsl_matrix_char_max
190 gsl_matrix_char_min
191 gsl_matrix_char_minmax
192 gsl_matrix_char_max_index
193 gsl_matrix_char_min_index
194 gsl_matrix_char_minmax_index
195 gsl_matrix_char_isnull
196 gsl_matrix_char_ispos
197 gsl_matrix_char_isneg
198 gsl_matrix_char_isnonneg
199 gsl_matrix_char_add
200 gsl_matrix_char_sub
201 gsl_matrix_char_mul_elements
202 gsl_matrix_char_div_elements
203 gsl_matrix_char_scale
204 gsl_matrix_char_add_constant
205 gsl_matrix_char_add_diagonal
208 double => [ qw/
209 gsl_matrix_alloc
210 gsl_matrix_calloc
211 gsl_matrix_alloc_from_block
212 gsl_matrix_alloc_from_matrix
213 gsl_vector_alloc_row_from_matrix
214 gsl_vector_alloc_col_from_matrix
215 gsl_matrix_free
216 gsl_matrix_submatrix
217 gsl_matrix_row
218 gsl_matrix_column
219 gsl_matrix_diagonal
220 gsl_matrix_subdiagonal
221 gsl_matrix_superdiagonal
222 gsl_matrix_subrow
223 gsl_matrix_subcolumn
224 gsl_matrix_view_array
225 gsl_matrix_view_array_with_tda
226 gsl_matrix_view_vector
227 gsl_matrix_view_vector_with_tda
228 gsl_matrix_const_submatrix
229 gsl_matrix_const_row
230 gsl_matrix_const_column
231 gsl_matrix_const_diagonal
232 gsl_matrix_const_subdiagonal
233 gsl_matrix_const_superdiagonal
234 gsl_matrix_const_subrow
235 gsl_matrix_const_subcolumn
236 gsl_matrix_const_view_array
237 gsl_matrix_const_view_array_with_tda
238 gsl_matrix_const_view_vector
239 gsl_matrix_const_view_vector_with_tda
240 gsl_matrix_get
241 gsl_matrix_set
242 gsl_matrix_ptr
243 gsl_matrix_const_ptr
244 gsl_matrix_set_zero
245 gsl_matrix_set_identity
246 gsl_matrix_set_all
247 gsl_matrix_fread
248 gsl_matrix_fwrite
249 gsl_matrix_fscanf
250 gsl_matrix_fprintf
251 gsl_matrix_memcpy
252 gsl_matrix_swap
253 gsl_matrix_swap_rows
254 gsl_matrix_swap_columns
255 gsl_matrix_swap_rowcol
256 gsl_matrix_transpose
257 gsl_matrix_transpose_memcpy
258 gsl_matrix_max
259 gsl_matrix_minmax
260 gsl_matrix_max_index
261 gsl_matrix_min_index
262 gsl_matrix_minmax_index
263 gsl_matrix_isnull
264 gsl_matrix_ispos
265 gsl_matrix_isneg
266 gsl_matrix_isnonneg
267 gsl_matrix_add
268 gsl_matrix_mul_elements
269 gsl_matrix_div_elements
270 gsl_matrix_scale
271 gsl_matrix_add_constant
272 gsl_matrix_add_diagonal
274 int => [ qw/
275 gsl_matrix_int_alloc
276 gsl_matrix_int_alloc_from_matrix
277 gsl_matrix_int_free
278 gsl_matrix_int_column
279 gsl_matrix_int_superdiagonal
280 gsl_matrix_int_view_array_with_tda
281 gsl_matrix_int_const_submatrix
282 gsl_matrix_int_const_diagonal
283 gsl_matrix_int_const_subrow
284 gsl_matrix_int_const_view_array_with_tda
285 gsl_matrix_int_get
286 gsl_matrix_int_ptr
287 gsl_matrix_int_set_zero
288 gsl_matrix_int_fread
289 gsl_matrix_int_fscanf
290 gsl_matrix_int_memcpy
291 gsl_matrix_int_swap_rows
292 gsl_matrix_int_transpose
293 gsl_matrix_int_max
294 gsl_matrix_int_max_index
295 gsl_matrix_int_minmax_index
296 gsl_matrix_int_ispos
297 gsl_matrix_int_add
298 gsl_matrix_int_mul_elements
299 gsl_matrix_int_add_constant
301 complex => [ qw/
302 gsl_matrix_complex_alloc
303 gsl_matrix_complex_calloc
304 gsl_matrix_complex_alloc_from_block
305 gsl_matrix_complex_alloc_from_matrix
306 gsl_vector_complex_alloc_row_from_matrix
307 gsl_vector_complex_alloc_col_from_matrix
308 gsl_matrix_complex_free
309 gsl_matrix_complex_submatrix
310 gsl_matrix_complex_row
311 gsl_matrix_complex_column
312 gsl_matrix_complex_diagonal
313 gsl_matrix_complex_subdiagonal
314 gsl_matrix_complex_superdiagonal
315 gsl_matrix_complex_subrow
316 gsl_matrix_complex_subcolumn
317 gsl_matrix_complex_view_array
318 gsl_matrix_complex_view_array_with_tda
319 gsl_matrix_complex_view_vector
320 gsl_matrix_complex_view_vector_with_tda
321 gsl_matrix_complex_const_submatrix
322 gsl_matrix_complex_const_row
323 gsl_matrix_complex_const_column
324 gsl_matrix_complex_const_diagonal
325 gsl_matrix_complex_const_subdiagonal
326 gsl_matrix_complex_const_superdiagonal
327 gsl_matrix_complex_const_subrow
328 gsl_matrix_complex_const_subcolumn
329 gsl_matrix_complex_const_view_array
330 gsl_matrix_complex_const_view_array_with_tda
331 gsl_matrix_complex_const_view_vector
332 gsl_matrix_complex_const_view_vector_with_tda
333 gsl_matrix_complex_get
334 gsl_matrix_complex_set
335 gsl_matrix_complex_ptr
336 gsl_matrix_complex_const_ptr
337 gsl_matrix_complex_set_zero gsl_matrix_complex_set_identity
338 gsl_matrix_complex_set_all
339 gsl_matrix_complex_fread
340 gsl_matrix_complex_fwrite
341 gsl_matrix_complex_fscanf
342 gsl_matrix_complex_fprintf
343 gsl_matrix_complex_memcpygsl_matrix_complex_swap
344 gsl_matrix_complex_swap_rows
345 gsl_matrix_complex_swap_columns
346 gsl_matrix_complex_swap_rowcol
347 gsl_matrix_complex_transpose
348 gsl_matrix_complex_transpose_memcpy
349 gsl_matrix_complex_isnull
350 gsl_matrix_complex_ispos
351 gsl_matrix_complex_isneg
352 gsl_matrix_complex_add
353 gsl_matrix_complex_sub
354 gsl_matrix_complex_mul_elements
355 gsl_matrix_complex_div_elements
356 gsl_matrix_complex_scale
357 gsl_matrix_complex_add_constant
358 gsl_matrix_complex_add_diagonal
359 gsl_matrix_complex_get_row
360 gsl_matrix_complex_get_col
361 gsl_matrix_complex_set_row
362 gsl_matrix_complex_set_col
363 /]);
365 =head1 NAME
367 Math::GSL::Matrix - Mathematical functions concerning Matrices
369 =head1 SYNOPSIS
371 use Math::GSL::Matrix qw/:all/;
372 my $matrix1 = Math::GSL::Matrix->new(5,5); # OO interface
373 my $matrix2 = gsl_matrix_alloc(5,5); # standard interface
376 =head1 Objected Oriented Interface to GSL Math::GSL::Matrix
378 =head2 Math::GSL::Matrix->new()
380 Creates a new Matrix of the given size.
382 my $matrix = Math::GSL::Matrix->new(10,10);
383 =cut
385 sub new
387 my ($class, $rows, $cols) = @_;
388 my $this = {};
389 my $matrix;
390 if ( defined $rows && defined $cols &&
391 $rows > 0 && $cols > 0 &&
392 (int $rows == $rows) && (int $cols == $cols)){
394 $matrix = gsl_matrix_alloc($rows,$cols);
395 } else {
396 croak( __PACKAGE__.'::new($x,$y) - $x and $y must be positive integers');
398 gsl_matrix_set_zero($matrix);
399 $this->{_matrix} = $matrix;
400 ($this->{_rows}, $this->{_cols}) = ($rows,$cols);
401 bless $this, $class;
403 =head2 raw()
405 Get the underlying GSL matrix object created by SWIG, useful for using gsl_matrix_* functions which do not have an OO counterpart.
407 my $matrix = Math::GSL::Matrix->new(3,3);
408 my $gsl_matrix = $matrix->raw;
409 my $stuff = gsl_matrix_get($gsl_matrix, 1, 2);
411 =cut
412 sub raw { (shift)->{_matrix} }
413 =head2 rows()
415 Returns the number of rows in the matrix.
417 my $rows = $matrix->rows;
418 =cut
420 sub rows { (shift)->{_rows} }
422 =head2 cols()
424 Returns the number of columns in the matrix.
426 my $cols = $matrix->cols;
427 =cut
429 sub cols { (shift)->{_cols} }
431 =head2 as_list()
433 Get the contents of a Math::GSL::Matrix object as a Perl list.
435 my $matrix = Math::GSL::Matrix->new(3,3);
437 my @matrix = $matrix->as_list;
438 =cut
441 sub as_list
443 my $self = shift;
444 my $line;
445 my @part;
446 my @total;
447 for($line=0; $line<$self->rows; $line++){
448 @part = map {
449 gsl_matrix_get($self->raw, $line, $_)
450 } (0 .. $self->cols-1 );
451 push(@total, @part);
453 return @total;
456 =head2 row()
458 Returns a row matrix of the row you enter.
460 my $matrix = Math::GSL::Matrix->new(3,3);
462 my $matrix_row = $matrix->row(0);
464 =cut
466 sub row
468 my ($self, $row) = @_;
469 croak (__PACKAGE__.'::$matrix->row($row) - invalid $row value')
470 unless (($row < $self->rows) and $row >= 0);
472 my $rowvec = Math::GSL::Vector->new($self->cols);
473 my $rowmat = Math::GSL::Matrix->new(1,$self->cols);
475 my $status = gsl_matrix_get_row($rowvec->raw, $self->raw, $row);
476 croak (__PACKAGE__.'::gsl_matrix_get_row - ' . gsl_strerror($status) )
477 unless ( $status == $GSL_SUCCESS );
479 $status = gsl_matrix_set_row($rowmat->raw, 0, $rowvec->raw);
481 croak (__PACKAGE__.'::gsl_matrix_set_row - ' . gsl_strerror($status) )
482 unless ( $status == $GSL_SUCCESS );
484 return $rowmat;
487 =head2 col()
489 Returns a col matrix of the column you enter.
491 my $matrix = Math::GSL::Matrix->new(3,3);
493 my $matrix_col = $matrix->col(0);
495 =cut
497 sub col
499 my ($self, $col) = @_;
500 croak (__PACKAGE__."::\$matrix->col(\$col) - $col not a valid column")
501 unless ($col < $self->cols and $col >= 0);
503 my $colvec = Math::GSL::Vector->new($self->cols);
504 my $colmat = Math::GSL::Matrix->new($self->rows, 1);
506 my $status = gsl_matrix_get_col($colvec->raw, $self->raw, $col);
507 # return $colvec;
508 $status = gsl_matrix_set_col($colmat->raw, 0, $colvec->raw);
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
821 The following functions are the same as the previous enonced ones except that they work with other data types than double.
823 =over 1
825 =item gsl_matrix_const_row gsl_matrix_const_column gsl_matrix_const_diagonal
827 =item gsl_matrix_const_subdiagonal gsl_matrix_const_superdiagonal
829 =item gsl_matrix_const_subrow gsl_matrix_const_subcolumn
831 =item gsl_matrix_const_view_array gsl_matrix_const_view_array_with_tda
833 =item gsl_matrix_const_view_vector gsl_matrix_const_view_vector_with_tda gsl_matrix_char_alloc gsl_matrix_char_calloc gsl_matrix_char_alloc_from_block
835 =item gsl_matrix_char_alloc_from_matrix gsl_vector_char_alloc_row_from_matrix gsl_vector_char_alloc_col_from_matrix
837 =item gsl_matrix_char_free gsl_matrix_char_submatrix gsl_matrix_char_row gsl_matrix_char_column
839 =item gsl_matrix_char_diagonal gsl_matrix_char_subdiagonal gsl_matrix_char_superdiagonal
841 =item gsl_matrix_char_subrow gsl_matrix_char_subcolumn gsl_matrix_char_view_array
843 =item gsl_matrix_char_view_array_with_tda gsl_matrix_char_view_vector gsl_matrix_char_view_vector_with_tda
845 =item gsl_matrix_char_const_submatrix gsl_matrix_char_const_row gsl_matrix_char_const_column
847 =item gsl_matrix_char_const_diagonal gsl_matrix_char_const_subdiagonal gsl_matrix_char_const_superdiagonal
849 =item gsl_matrix_char_const_subrow gsl_matrix_char_const_subcolumn gsl_matrix_char_const_view_array
851 =item gsl_matrix_char_const_view_array_with_tda gsl_matrix_char_const_view_vector gsl_matrix_char_const_view_vector_with_tda
853 =item gsl_matrix_char_get gsl_matrix_char_set gsl_matrix_char_ptr gsl_matrix_char_const_ptr
855 =item gsl_matrix_char_set_zero gsl_matrix_char_set_identity
857 =item gsl_matrix_char_set_all gsl_matrix_char_fread
859 =item gsl_matrix_char_fwrite gsl_matrix_char_fscanf gsl_matrix_char_fprintf
861 =item gsl_matrix_char_memcpy gsl_matrix_char_swap
863 =item gsl_matrix_char_swap_rows gsl_matrix_char_swap_columns
865 =item gsl_matrix_char_swap_rowcol gsl_matrix_char_transpose gsl_matrix_char_transpose_memcpy
867 =item gsl_matrix_char_max gsl_matrix_char_min
869 =item gsl_matrix_char_minmax gsl_matrix_char_max_index
871 =item gsl_matrix_char_min_index gsl_matrix_char_minmax_index
873 =item gsl_matrix_char_isnull gsl_matrix_char_ispos gsl_matrix_char_isneg
875 =item gsl_matrix_char_isnonneg gsl_matrix_char_add
877 =item gsl_matrix_char_sub gsl_matrix_char_mul_elements gsl_matrix_char_div_elements
879 =item gsl_matrix_char_scale gsl_matrix_char_add_constant gsl_matrix_char_add_diagonal
881 =item gsl_matrix_int_alloc gsl_matrix_int_calloc gsl_matrix_int_alloc_from_block
883 =item gsl_matrix_int_alloc_from_matrix gsl_vector_int_alloc_row_from_matrix gsl_vector_int_alloc_col_from_matrix
885 =item gsl_matrix_int_free gsl_matrix_int_submatrix gsl_matrix_int_row
887 =item gsl_matrix_int_column gsl_matrix_int_diagonal gsl_matrix_int_subdiagonal gsl_matrix_int_superdiagonal gsl_matrix_int_subrow gsl_matrix_int_subcolumn gsl_matrix_int_view_array gsl_matrix_int_view_array_with_tda gsl_matrix_int_view_vector gsl_matrix_int_view_vector_with_tda
889 =item gsl_matrix_int_const_submatrix gsl_matrix_int_const_row gsl_matrix_int_const_column
892 =item gsl_matrix_int_const_diagonal gsl_matrix_int_const_subdiagonal gsl_matrix_int_const_superdiagonal
894 =item gsl_matrix_int_const_subrow gsl_matrix_int_const_subcolumn gsl_matrix_int_const_view_array
896 =item gsl_matrix_int_const_view_array_with_tda gsl_matrix_int_const_view_vector gsl_matrix_int_const_view_vector_with_tda
898 =item gsl_matrix_int_get gsl_matrix_int_set
900 =item gsl_matrix_int_ptr gsl_matrix_int_const_ptr
902 =item gsl_matrix_int_set_zero gsl_matrix_int_set_identity gsl_matrix_int_set_all
904 =item gsl_matrix_int_fread gsl_matrix_int_fwrite
906 =item gsl_matrix_int_fscanf gsl_matrix_int_fprintf
908 =item gsl_matrix_int_memcpy gsl_matrix_int_swap
910 =item gsl_matrix_int_swap_rows gsl_matrix_int_swap_columns gsl_matrix_int_swap_rowcol
912 =item gsl_matrix_int_transpose gsl_matrix_int_transpose_memcpy
914 =item gsl_matrix_int_max gsl_matrix_int_min gsl_matrix_int_minmax
916 =item gsl_matrix_int_max_index gsl_matrix_int_min_index
918 =item gsl_matrix_int_minmax_index gsl_matrix_int_isnull
920 =item gsl_matrix_int_ispos gsl_matrix_int_isneg gsl_matrix_int_isnonneg
922 =item gsl_matrix_int_add gsl_matrix_int_sub
924 =item gsl_matrix_int_mul_elements gsl_matrix_int_div_elements gsl_matrix_int_scale
926 =item gsl_matrix_int_add_constant gsl_matrix_int_add_diagonal
928 =back
930 You have to add the functions you want to use inside the qw /put_funtion_here /.
931 You can also write use Math::GSL::Matrix qw/:all/ to use all avaible functions of the module.
932 Other tags are also avaible, here is a complete list of all tags for this module :
934 =over 1
936 =item C<all>
938 =item C<int>
940 =item C<double>
942 =item C<char>
944 =item C<complex>
946 =back
948 For more informations on the functions, we refer you to the GSL offcial documentation: http://www.gnu.org/software/gsl/manual/html_node/
949 Tip : search on google: site:http://www.gnu.org/software/gsl/manual/html_node/name_of_the_function_you_want
951 =head1 EXAMPLES
953 Most of the examples from this section are perl versions of the examples there : http://www.gnu.org/software/gsl/manual/html_node/Example-programs-for-matrices.html
955 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.
957 use Math::GSL::Matrix qw/:all/;
958 $m = gsl_matrix_alloc (10,3);
959 for($i = 0; $i< 10; $i++){
960 for($j = 0; $j<3; $j++){
961 gsl_matrix_set($m, $i, $j, 0.23 + 100*$i + $j);
965 for($i = 0; $i< 100; $i++){ # OUT OF RANGE ERROR
966 for($j=0; $j<3; $j++){
967 print "m($i, $j) = " . gsl_matrix_get ($m, $i, $j) . "\n";
970 gsl_matrix_free ($m);
973 use Math::GSL::Matrix qw/:all/;
975 $m = gsl_matrix_alloc (100, 100);
976 $a = gsl_matrix_alloc (100, 100);
978 for($i = 0; $i < 100; $i++){
979 for($j = 0; $j < 100; $j++){
980 gsl_matrix_set ($m, $i, $j, 0.23 + $i + $j);
984 The next program shows how to write a matrix to a file.
986 $out = fopen("test.dat", "wb");
987 gsl_matrix_fwrite ($out, $m);
988 fclose ($out);
990 $in = fopen("test.dat", "rb");
991 gsl_matrix_fread ($in, $a);
992 fclose($in);
994 $k=0;
995 for($i = 0; $i < 100; $i++){
996 for($j = 0; $j < 100; $j++){
997 $mij = gsl_matrix_get ($m, $i, $j);
998 $aij = gsl_matrix_get ($a, $i, $j);
999 if ($mij != $aij){ $k++ };
1003 gsl_matrix_free($m);
1004 gsl_matrix_free($a);
1006 print "differences = $k (should be zero)\n";
1008 =head1 AUTHOR
1010 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
1012 =head1 COPYRIGHT AND LICENSE
1014 Copyright (C) 2008 Jonathan Leto and Thierry Moisan
1016 This program is free software; you can redistribute it and/or modify it
1017 under the same terms as Perl itself.
1019 =cut