Adding an example to the Fit documentation
[Math-GSL.git] / Matrix.i
blob15b8cc5b42aebb9b4ed4ddb16795e6fd94a368b4
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 gsl_vector *OUTPUT { gsl_vector *V };
8 %apply double *OUTPUT { double * min_out, double * max_out };
10 FILE * fopen(char *, char *);
11 int fclose(FILE *);
14 #include "gsl/gsl_matrix.h"
15 #include "gsl/gsl_complex.h"
16 #include "gsl/gsl_vector_double.h"
17 #include "gsl/gsl_matrix_double.h"
18 #include "gsl/gsl_matrix_int.h"
19 #include "gsl/gsl_matrix_complex_double.h"
20 #include "gsl/gsl_matrix_char.h"
23 %include "gsl/gsl_matrix.h"
24 %include "gsl/gsl_complex.h"
25 %include "gsl/gsl_vector_double.h"
26 %include "gsl/gsl_matrix_double.h"
27 %include "gsl/gsl_matrix_int.h"
28 %include "gsl/gsl_matrix_complex_double.h"
29 %include "gsl/gsl_matrix_char.h"
32 %perlcode %{
34 no warnings 'redefine';
35 use Carp qw/croak/;
36 use Math::GSL qw/:all/;
37 use Math::GSL::Errno qw/:all/;
39 @EXPORT_OK = qw/fopen fclose
40 gsl_matrix_alloc gsl_matrix_calloc gsl_matrix_alloc_from_block
41 gsl_matrix_alloc_from_matrix gsl_vector_alloc_row_from_matrix
42 gsl_vector_alloc_col_from_matrix gsl_matrix_free gsl_matrix_submatrix
43 gsl_matrix_row gsl_matrix_column gsl_matrix_diagonal
44 gsl_matrix_subdiagonal gsl_matrix_superdiagonal gsl_matrix_subrow
45 gsl_matrix_subcolumn gsl_matrix_view_array
46 gsl_matrix_view_array_with_tda gsl_matrix_view_vector
47 gsl_matrix_view_vector_with_tda gsl_matrix_const_submatrix
48 gsl_matrix_const_row gsl_matrix_const_column gsl_matrix_const_diagonal
49 gsl_matrix_const_subdiagonal gsl_matrix_const_superdiagonal
50 gsl_matrix_const_subrow gsl_matrix_const_subcolumn
51 gsl_matrix_const_view_array gsl_matrix_const_view_array_with_tda
52 gsl_matrix_const_view_vector gsl_matrix_const_view_vector_with_tda
53 gsl_matrix_get gsl_matrix_set gsl_matrix_ptr gsl_matrix_const_ptr
54 gsl_matrix_set_zero gsl_matrix_set_identity gsl_matrix_set_all
55 gsl_matrix_fread gsl_matrix_fwrite gsl_matrix_fscanf gsl_matrix_fprintf
56 gsl_matrix_memcpy gsl_matrix_swap gsl_matrix_swap_rows
57 gsl_matrix_swap_columns gsl_matrix_swap_rowcol gsl_matrix_transpose
58 gsl_matrix_transpose_memcpy gsl_matrix_max gsl_matrix_min
59 gsl_matrix_minmax gsl_matrix_max_index gsl_matrix_min_index
60 gsl_matrix_minmax_index gsl_matrix_isnull gsl_matrix_ispos
61 gsl_matrix_isneg gsl_matrix_isnonneg gsl_matrix_add gsl_matrix_sub
62 gsl_matrix_mul_elements gsl_matrix_div_elements gsl_matrix_scale
63 gsl_matrix_add_constant gsl_matrix_add_diagonal
64 gsl_matrix_char_alloc gsl_matrix_char_calloc gsl_matrix_char_alloc_from_block
65 gsl_matrix_char_alloc_from_matrix gsl_vector_char_alloc_row_from_matrix gsl_vector_char_alloc_col_from_matrix
66 gsl_matrix_char_free gsl_matrix_char_submatrix
67 gsl_matrix_char_row gsl_matrix_char_column
68 gsl_matrix_char_diagonal gsl_matrix_char_subdiagonal gsl_matrix_char_superdiagonal
69 gsl_matrix_char_subrow gsl_matrix_char_subcolumn gsl_matrix_char_view_array
70 gsl_matrix_char_view_array_with_tda gsl_matrix_char_view_vector gsl_matrix_char_view_vector_with_tda
71 gsl_matrix_char_const_submatrix gsl_matrix_char_const_row gsl_matrix_char_const_column
72 gsl_matrix_char_const_diagonal gsl_matrix_char_const_subdiagonal gsl_matrix_char_const_superdiagonal
73 gsl_matrix_char_const_subrow gsl_matrix_char_const_subcolumn gsl_matrix_char_const_view_array
74 gsl_matrix_char_const_view_array_with_tda gsl_matrix_char_const_view_vector gsl_matrix_char_const_view_vector_with_tda
75 gsl_matrix_char_get gsl_matrix_char_set gsl_matrix_char_ptr gsl_matrix_char_const_ptr
76 gsl_matrix_char_set_zero gsl_matrix_char_set_identity
77 gsl_matrix_char_set_all gsl_matrix_char_fread
78 gsl_matrix_char_fwrite gsl_matrix_char_fscanf gsl_matrix_char_fprintf
79 gsl_matrix_char_memcpy gsl_matrix_char_swap
80 gsl_matrix_char_swap_rows gsl_matrix_char_swap_columns
81 gsl_matrix_char_swap_rowcol gsl_matrix_char_transpose gsl_matrix_char_transpose_memcpy
82 gsl_matrix_char_max gsl_matrix_char_min
83 gsl_matrix_char_minmax gsl_matrix_char_max_index
84 gsl_matrix_char_min_index gsl_matrix_char_minmax_index
85 gsl_matrix_char_isnull gsl_matrix_char_ispos gsl_matrix_char_isneg
86 gsl_matrix_char_isnonneg gsl_matrix_char_add
87 gsl_matrix_char_sub gsl_matrix_char_mul_elements gsl_matrix_char_div_elements
88 gsl_matrix_char_scale gsl_matrix_char_add_constant gsl_matrix_char_add_diagonal
89 gsl_matrix_int_alloc gsl_matrix_int_calloc gsl_matrix_int_alloc_from_block
90 gsl_matrix_int_alloc_from_matrix gsl_vector_int_alloc_row_from_matrix gsl_vector_int_alloc_col_from_matrix
91 gsl_matrix_int_free gsl_matrix_int_submatrix gsl_matrix_int_row
92 gsl_matrix_int_column gsl_matrix_int_diagonal gsl_matrix_int_subdiagonal
93 gsl_matrix_int_superdiagonal gsl_matrix_int_subrow gsl_matrix_int_subcolumn gsl_matrix_int_view_array
94 gsl_matrix_int_view_array_with_tda gsl_matrix_int_view_vector gsl_matrix_int_view_vector_with_tda
95 gsl_matrix_int_const_submatrix gsl_matrix_int_const_row gsl_matrix_int_const_column
96 gsl_matrix_int_const_diagonal gsl_matrix_int_const_subdiagonal gsl_matrix_int_const_superdiagonal
97 gsl_matrix_int_const_subrow gsl_matrix_int_const_subcolumn gsl_matrix_int_const_view_array
98 gsl_matrix_int_const_view_array_with_tda gsl_matrix_int_const_view_vector gsl_matrix_int_const_view_vector_with_tda
99 gsl_matrix_int_get gsl_matrix_int_set
100 gsl_matrix_int_ptr gsl_matrix_int_const_ptr
101 gsl_matrix_int_set_zero gsl_matrix_int_set_identity gsl_matrix_int_set_all
102 gsl_matrix_int_fread gsl_matrix_int_fwrite
103 gsl_matrix_int_fscanf gsl_matrix_int_fprintf
104 gsl_matrix_int_memcpy gsl_matrix_int_swap
105 gsl_matrix_int_swap_rows gsl_matrix_int_swap_columns gsl_matrix_int_swap_rowcol
106 gsl_matrix_int_transpose gsl_matrix_int_transpose_memcpy
107 gsl_matrix_int_max gsl_matrix_int_min gsl_matrix_int_minmax
108 gsl_matrix_int_max_index gsl_matrix_int_min_index
109 gsl_matrix_int_minmax_index gsl_matrix_int_isnull
110 gsl_matrix_int_ispos gsl_matrix_int_isneg gsl_matrix_int_isnonneg
111 gsl_matrix_int_add gsl_matrix_int_sub
112 gsl_matrix_int_mul_elements gsl_matrix_int_div_elements gsl_matrix_int_scale
113 gsl_matrix_int_add_constant gsl_matrix_int_add_diagonal
114 gsl_matrix_get_row gsl_matrix_get_col gsl_matrix_set_row gsl_matrix_set_col
115 gsl_matrix_complex_alloc gsl_matrix_complex_calloc gsl_matrix_complex_alloc_from_block
116 gsl_matrix_complex_alloc_from_matrix gsl_vector_complex_alloc_row_from_matrix gsl_vector_complex_alloc_col_from_matrix
117 gsl_matrix_complex_free gsl_matrix_complex_submatrix gsl_matrix_complex_row
118 gsl_matrix_complex_column gsl_matrix_complex_diagonal gsl_matrix_complex_subdiagonal
119 gsl_matrix_complex_superdiagonal gsl_matrix_complex_subrow gsl_matrix_complex_subcolumn
120 gsl_matrix_complex_view_array gsl_matrix_complex_view_array_with_tda gsl_matrix_complex_view_vector
121 gsl_matrix_complex_view_vector_with_tda gsl_matrix_complex_const_submatrix gsl_matrix_complex_const_row
122 gsl_matrix_complex_const_column gsl_matrix_complex_const_diagonal gsl_matrix_complex_const_subdiagonal
123 gsl_matrix_complex_const_superdiagonal gsl_matrix_complex_const_subrow gsl_matrix_complex_const_subcolumn
124 gsl_matrix_complex_const_view_array gsl_matrix_complex_const_view_array_with_tda gsl_matrix_complex_const_view_vector
125 gsl_matrix_complex_const_view_vector_with_tda gsl_matrix_complex_get gsl_matrix_complex_set
126 gsl_matrix_complex_ptr gsl_matrix_complex_const_ptr gsl_matrix_complex_set_zero
127 gsl_matrix_complex_set_identity gsl_matrix_complex_set_all gsl_matrix_complex_fread
128 gsl_matrix_complex_fwrite gsl_matrix_complex_fscanf gsl_matrix_complex_fprintf
129 gsl_matrix_complex_memcpy gsl_matrix_complex_swap gsl_matrix_complex_swap_rows
130 gsl_matrix_complex_swap_columns gsl_matrix_complex_swap_rowcol gsl_matrix_complex_transpose
131 gsl_matrix_complex_transpose_memcpy gsl_matrix_complex_isnull gsl_matrix_complex_ispos
132 gsl_matrix_complex_isneg gsl_matrix_complex_add gsl_matrix_complex_sub
133 gsl_matrix_complex_mul_elements gsl_matrix_complex_div_elements gsl_matrix_complex_scale
134 gsl_matrix_complex_add_constant gsl_matrix_complex_add_diagonal gsl_matrix_complex_get_row
135 gsl_matrix_complex_get_col gsl_matrix_complex_set_row gsl_matrix_complex_set_col /;
139 %EXPORT_TAGS = ( all => [ @EXPORT_OK ],
140 char => [ qw/
141 gsl_matrix_char_alloc
142 gsl_matrix_char_calloc
143 gsl_matrix_char_alloc_from_block
144 gsl_matrix_char_alloc_from_matrix
145 gsl_vector_char_alloc_row_from_matrix
146 gsl_vector_char_alloc_col_from_matrix
147 gsl_matrix_char_free
148 gsl_matrix_char_submatrix
149 gsl_matrix_char_row
150 gsl_matrix_char_column
151 gsl_matrix_char_diagonal
152 gsl_matrix_char_subdiagonal
153 gsl_matrix_char_superdiagonal
154 gsl_matrix_char_subrow
155 gsl_matrix_char_subcolumn
156 gsl_matrix_char_view_array
157 gsl_matrix_char_view_array_with_tda
158 gsl_matrix_char_view_vector
159 gsl_matrix_char_view_vector_with_tda
160 gsl_matrix_char_const_submatrix
161 gsl_matrix_char_const_row
162 gsl_matrix_char_const_column
163 gsl_matrix_char_const_diagonal
164 gsl_matrix_char_const_subdiagonal
165 gsl_matrix_char_const_superdiagonal
166 gsl_matrix_char_const_subrow
167 gsl_matrix_char_const_subcolumn
168 gsl_matrix_char_const_view_array
169 gsl_matrix_char_const_view_array_with_tda
170 gsl_matrix_char_const_view_vector
171 gsl_matrix_char_const_view_vector_with_tda
172 gsl_matrix_char_get
173 gsl_matrix_char_set
174 gsl_matrix_char_ptr
175 gsl_matrix_char_const_ptr
176 gsl_matrix_char_set_zero
177 gsl_matrix_char_set_identity
178 gsl_matrix_char_set_all
179 gsl_matrix_char_fread
180 gsl_matrix_char_fwrite
181 gsl_matrix_char_fscanf
182 gsl_matrix_char_fprintf
183 gsl_matrix_char_memcpy
184 gsl_matrix_char_swap
185 gsl_matrix_char_swap_rows
186 gsl_matrix_char_swap_columns
187 gsl_matrix_char_swap_rowcol
188 gsl_matrix_char_transpose
189 gsl_matrix_char_transpose_memcpy
190 gsl_matrix_char_max
191 gsl_matrix_char_min
192 gsl_matrix_char_minmax
193 gsl_matrix_char_max_index
194 gsl_matrix_char_min_index
195 gsl_matrix_char_minmax_index
196 gsl_matrix_char_isnull
197 gsl_matrix_char_ispos
198 gsl_matrix_char_isneg
199 gsl_matrix_char_isnonneg
200 gsl_matrix_char_add
201 gsl_matrix_char_sub
202 gsl_matrix_char_mul_elements
203 gsl_matrix_char_div_elements
204 gsl_matrix_char_scale
205 gsl_matrix_char_add_constant
206 gsl_matrix_char_add_diagonal
209 double => [ qw/
210 gsl_matrix_alloc
211 gsl_matrix_calloc
212 gsl_matrix_alloc_from_block
213 gsl_matrix_alloc_from_matrix
214 gsl_vector_alloc_row_from_matrix
215 gsl_vector_alloc_col_from_matrix
216 gsl_matrix_free
217 gsl_matrix_submatrix
218 gsl_matrix_row
219 gsl_matrix_column
220 gsl_matrix_diagonal
221 gsl_matrix_subdiagonal
222 gsl_matrix_superdiagonal
223 gsl_matrix_subrow
224 gsl_matrix_subcolumn
225 gsl_matrix_view_array
226 gsl_matrix_view_array_with_tda
227 gsl_matrix_view_vector
228 gsl_matrix_view_vector_with_tda
229 gsl_matrix_const_submatrix
230 gsl_matrix_const_row
231 gsl_matrix_const_column
232 gsl_matrix_const_diagonal
233 gsl_matrix_const_subdiagonal
234 gsl_matrix_const_superdiagonal
235 gsl_matrix_const_subrow
236 gsl_matrix_const_subcolumn
237 gsl_matrix_const_view_array
238 gsl_matrix_const_view_array_with_tda
239 gsl_matrix_const_view_vector
240 gsl_matrix_const_view_vector_with_tda
241 gsl_matrix_get
242 gsl_matrix_set
243 gsl_matrix_ptr
244 gsl_matrix_const_ptr
245 gsl_matrix_set_zero
246 gsl_matrix_set_identity
247 gsl_matrix_set_all
248 gsl_matrix_fread
249 gsl_matrix_fwrite
250 gsl_matrix_fscanf
251 gsl_matrix_fprintf
252 gsl_matrix_memcpy
253 gsl_matrix_swap
254 gsl_matrix_swap_rows
255 gsl_matrix_swap_columns
256 gsl_matrix_swap_rowcol
257 gsl_matrix_transpose
258 gsl_matrix_transpose_memcpy
259 gsl_matrix_max
260 gsl_matrix_minmax
261 gsl_matrix_max_index
262 gsl_matrix_min_index
263 gsl_matrix_minmax_index
264 gsl_matrix_isnull
265 gsl_matrix_ispos
266 gsl_matrix_isneg
267 gsl_matrix_isnonneg
268 gsl_matrix_add
269 gsl_matrix_mul_elements
270 gsl_matrix_div_elements
271 gsl_matrix_scale
272 gsl_matrix_add_constant
273 gsl_matrix_add_diagonal
275 int => [ qw/
276 gsl_matrix_int_alloc
277 gsl_matrix_int_alloc_from_matrix
278 gsl_matrix_int_free
279 gsl_matrix_int_column
280 gsl_matrix_int_superdiagonal
281 gsl_matrix_int_view_array_with_tda
282 gsl_matrix_int_const_submatrix
283 gsl_matrix_int_const_diagonal
284 gsl_matrix_int_const_subrow
285 gsl_matrix_int_const_view_array_with_tda
286 gsl_matrix_int_get
287 gsl_matrix_int_ptr
288 gsl_matrix_int_set_zero
289 gsl_matrix_int_fread
290 gsl_matrix_int_fscanf
291 gsl_matrix_int_memcpy
292 gsl_matrix_int_swap_rows
293 gsl_matrix_int_transpose
294 gsl_matrix_int_max
295 gsl_matrix_int_max_index
296 gsl_matrix_int_minmax_index
297 gsl_matrix_int_ispos
298 gsl_matrix_int_add
299 gsl_matrix_int_mul_elements
300 gsl_matrix_int_add_constant
302 complex => [ qw/
303 gsl_matrix_complex_alloc
304 gsl_matrix_complex_calloc
305 gsl_matrix_complex_alloc_from_block
306 gsl_matrix_complex_alloc_from_matrix
307 gsl_vector_complex_alloc_row_from_matrix
308 gsl_vector_complex_alloc_col_from_matrix
309 gsl_matrix_complex_free
310 gsl_matrix_complex_submatrix
311 gsl_matrix_complex_row
312 gsl_matrix_complex_column
313 gsl_matrix_complex_diagonal
314 gsl_matrix_complex_subdiagonal
315 gsl_matrix_complex_superdiagonal
316 gsl_matrix_complex_subrow
317 gsl_matrix_complex_subcolumn
318 gsl_matrix_complex_view_array
319 gsl_matrix_complex_view_array_with_tda
320 gsl_matrix_complex_view_vector
321 gsl_matrix_complex_view_vector_with_tda
322 gsl_matrix_complex_const_submatrix
323 gsl_matrix_complex_const_row
324 gsl_matrix_complex_const_column
325 gsl_matrix_complex_const_diagonal
326 gsl_matrix_complex_const_subdiagonal
327 gsl_matrix_complex_const_superdiagonal
328 gsl_matrix_complex_const_subrow
329 gsl_matrix_complex_const_subcolumn
330 gsl_matrix_complex_const_view_array
331 gsl_matrix_complex_const_view_array_with_tda
332 gsl_matrix_complex_const_view_vector
333 gsl_matrix_complex_const_view_vector_with_tda
334 gsl_matrix_complex_get
335 gsl_matrix_complex_set
336 gsl_matrix_complex_ptr
337 gsl_matrix_complex_const_ptr
338 gsl_matrix_complex_set_zero gsl_matrix_complex_set_identity
339 gsl_matrix_complex_set_all
340 gsl_matrix_complex_fread
341 gsl_matrix_complex_fwrite
342 gsl_matrix_complex_fscanf
343 gsl_matrix_complex_fprintf
344 gsl_matrix_complex_memcpygsl_matrix_complex_swap
345 gsl_matrix_complex_swap_rows
346 gsl_matrix_complex_swap_columns
347 gsl_matrix_complex_swap_rowcol
348 gsl_matrix_complex_transpose
349 gsl_matrix_complex_transpose_memcpy
350 gsl_matrix_complex_isnull
351 gsl_matrix_complex_ispos
352 gsl_matrix_complex_isneg
353 gsl_matrix_complex_add
354 gsl_matrix_complex_sub
355 gsl_matrix_complex_mul_elements
356 gsl_matrix_complex_div_elements
357 gsl_matrix_complex_scale
358 gsl_matrix_complex_add_constant
359 gsl_matrix_complex_add_diagonal
360 gsl_matrix_complex_get_row
361 gsl_matrix_complex_get_col
362 gsl_matrix_complex_set_row
363 gsl_matrix_complex_set_col
364 /]);
366 =head1 NAME
368 Math::GSL::Matrix - Mathematical functions concerning Matrices
370 =head1 SYNOPSIS
372 use Math::GSL::Matrix qw/:all/;
373 my $matrix1 = Math::GSL::Matrix->new(5,5); # OO interface
374 my $matrix2 = gsl_matrix_alloc(5,5); # standard interface
377 =head1 Objected Oriented Interface to GSL Math::GSL::Matrix
379 =head2 Math::GSL::Matrix->new()
381 Creates a new Matrix of the given size.
383 my $matrix = Math::GSL::Matrix->new(10,10);
384 =cut
386 sub new
388 my ($class, $rows, $cols) = @_;
389 my $this = {};
390 my $matrix;
391 if ( defined $rows && defined $cols &&
392 $rows > 0 && $cols > 0 &&
393 (int $rows == $rows) && (int $cols == $cols)){
395 $matrix = gsl_matrix_alloc($rows,$cols);
396 } else {
397 croak( __PACKAGE__.'::new($x,$y) - $x and $y must be positive integers');
399 gsl_matrix_set_zero($matrix);
400 $this->{_matrix} = $matrix;
401 ($this->{_rows}, $this->{_cols}) = ($rows,$cols);
402 bless $this, $class;
404 =head2 raw()
406 Get the underlying GSL matrix object created by SWIG, useful for using gsl_matrix_* functions which do not have an OO counterpart.
408 my $matrix = Math::GSL::Matrix->new(3,3);
409 my $gsl_matrix = $matrix->raw;
410 my $stuff = gsl_matrix_get($gsl_matrix, 1, 2);
412 =cut
413 sub raw { (shift)->{_matrix} }
414 =head2 rows()
416 Returns the number of rows in the matrix.
418 my $rows = $matrix->rows;
419 =cut
421 sub rows { (shift)->{_rows} }
423 =head2 cols()
425 Returns the number of columns in the matrix.
427 my $cols = $matrix->cols;
428 =cut
430 sub cols { (shift)->{_cols} }
432 =head2 as_list()
434 Get the contents of a Math::GSL::Matrix object as a Perl list.
436 my $matrix = Math::GSL::Matrix->new(3,3);
438 my @matrix = $matrix->as_list;
439 =cut
442 sub as_list
444 my $self = shift;
445 my $line;
446 my @part;
447 my @total;
448 for($line=0; $line<$self->rows; $line++){
449 @part = map {
450 gsl_matrix_get($self->raw, $line, $_)
451 } (0 .. $self->cols-1 );
452 push(@total, @part);
454 return @total;
457 =head2 row()
459 Returns a row matrix of the row you enter.
461 my $matrix = Math::GSL::Matrix->new(3,3);
463 my $matrix_row = $matrix->row(0);
465 =cut
467 sub row
469 my ($self, $row) = @_;
470 croak (__PACKAGE__.'::$matrix->row($row) - invalid $row value')
471 unless (($row < $self->rows) and $row >= 0);
473 my $rowvec = Math::GSL::Vector->new($self->cols);
474 my $rowmat = Math::GSL::Matrix->new(1,$self->cols);
476 my $status = gsl_matrix_get_row($rowvec->raw, $self->raw, $row);
477 croak (__PACKAGE__.'::gsl_matrix_get_row - ' . gsl_strerror($status) )
478 unless ( $status == $GSL_SUCCESS );
480 $status = gsl_matrix_set_row($rowmat->raw, 0, $rowvec->raw);
482 croak (__PACKAGE__.'::gsl_matrix_set_row - ' . gsl_strerror($status) )
483 unless ( $status == $GSL_SUCCESS );
485 return $rowmat;
488 =head2 col()
490 Returns a col matrix of the column you enter.
492 my $matrix = Math::GSL::Matrix->new(3,3);
494 my $matrix_col = $matrix->col(0);
496 =cut
498 sub col
500 my ($self, $col) = @_;
501 croak (__PACKAGE__."::\$matrix->col(\$col) - $col not a valid column")
502 unless ($col < $self->cols and $col >= 0);
504 my $colvec = Math::GSL::Vector->new($self->cols);
505 my $colmat = Math::GSL::Matrix->new($self->rows, 1);
507 my $status = gsl_matrix_get_col($colvec->raw, $self->raw, $col);
508 # return $colvec;
509 $status = gsl_matrix_set_col($colmat->raw, 0, $colvec->raw);
510 return $colmat;
513 =head2 set_row()
515 Sets a the values of a row with the elements of an array.
517 my $matrix = Math::GSL::Matrix->new(3,3);
518 $matrix->set_row(0, [8, 6, 2]);
520 You can also set multiple rows at once with chained calls:
521 my $matrix = Math::GSL::Matrix->new(3,3);
522 $matrix->set_row(0, [8, 6, 2]);
523 ->set_row(1, [2, 4, 1]);
526 =cut
528 sub set_row {
529 my ($self, $row, $values) = @_;
530 my $length = $#$values;
531 die __PACKAGE__.'::new($x, $values) - $values must be a nonempty array reference' if $length == -1;
532 die __PACKAGE__.'::set_row($x, $values) - $x must be a valid row number' if ($row < 0 || $row >= $self->rows);
533 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);
534 map { gsl_matrix_set($self->raw, $row, $_, $values->[$_]) } (0..$length);
535 return $self;
538 =head2 set_col()
540 Sets a the values of a column with the elements of an array.
542 my $matrix = Math::GSL::Matrix->new(3,3);
543 $matrix->set_col(0, [8, 6, 2]);
545 You can also set multiple columns at once with chained calls:
546 my $matrix = Math::GSL::Matrix->new(3,3);
547 $matrix->set_col(0, [8, 6, 2]);
548 ->set_col(1, [2, 4, 1]);
551 =cut
553 sub set_col {
554 my ($self, $col, $values) = @_;
555 my $length = $#$values;
556 die __PACKAGE__.'::new($x, $values) - $values must be a nonempty array reference' if $length == -1;
557 die __PACKAGE__.'::set_row($x, $values) - $x must be a valid column number' if ($col < 0 || $col >= $self->cols);
558 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);
559 map { gsl_matrix_set($self->raw, $_, $col, $values->[$_]) } (0..$length);
560 return $self;
563 =head1 DESCRIPTION
565 Here is a list of all the functions included in this module :
567 =over 1
569 =item C<gsl_matrix_alloc($i, $j)> - Return a gsl_matrix of $i rows and $j columns
571 =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
573 =item C<gsl_matrix_alloc_from_block> -
575 =item C<gsl_matrix_free> -
577 =item C<gsl_matrix_alloc_from_matrix > -
579 =item C<gsl_vector_alloc_row_from_matrix> -
581 =item C<gsl_vector_alloc_col_from_matrix > -
583 =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.
585 =item C<gsl_matrix_row($m , $i)> - Return a vector view of the $i-th row of the matrix $m
587 =item C<gsl_matrix_column($m, $j)> - Return a vector view of the $j-th column of the matrix $m
589 =item C<gsl_matrix_diagonal($m)> - Return a vector view of the diagonal of the vector. The matrix doesn't have to be square.
591 =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.
593 =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.
595 =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.
597 =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.
599 =item C<gsl_matrix_view_array> -
601 =item C<gsl_matrix_view_array_with_tda > -
603 =item C<gsl_matrix_view_vector> -
605 =item C<gsl_matrix_view_vector_with_tda> -
607 =item C<gsl_matrix_const_submatrix> -
609 =item C<gsl_matrix_get($m, $i, $j)> - Return the (i,j)-th element of the matrix $m
611 =item C<gsl_matrix_set($m, $i, $j, $x)> - Set the value of the (i,j)-th element of the matrix $m to $x
613 =item C<gsl_matrix_ptr> -
615 =item C<gsl_matrix_const_ptr> -
617 =item C<gsl_matrix_set_zero($m)> - Set all the elements of the matrix $m to zero
619 =item C<gsl_matrix_set_identity($m)> - Set the elements of the matrix $m to the corresponding elements of the identity matrix
621 =item C<gsl_matrix_set_all($m, $x)> - Set all the elements of the matrix $m to the value $x
623 =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
625 =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
627 =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
629 =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
631 =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.
633 =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.
635 =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.
637 =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.
639 =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.
641 =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.
643 =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.
645 =item C<gsl_matrix_max($m)> - Return the maximum value in the matrix $m
647 =item C<gsl_matrix_min($m)> - Return the minimum value in the matrix $m
649 =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.
651 =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.
653 =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.
655 =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
657 =item C<gsl_matrix_isnull($m)> - Return 1 if all the elements of the matrix $m are zero, 0 otherwise
659 =item C<gsl_matrix_ispos($m)> - Return 1 if all the elements of the matrix $m are strictly positve, 0 otherwise
661 =item C<gsl_matrix_isneg($m)> - Return 1 if all the elements of the matrix $m are strictly negative, 0 otherwise
663 =item C<gsl_matrix_isnonneg($m)> - Return 1 if all the elements of the matrix $m are non-negatuive, 0 otherwise
665 =item C<gsl_matrix_add($a, $b)> - Add the elements of matrix $b to the elements of matrix $a
667 =item C<gsl_matrix_sub($a, $b)> - Subtract the elements of matrix $b from the elements of matrix $a
669 =item C<gsl_matrix_mul_elements($a, $b)> - Multiplie the elements of matrix $a by the elements of matrix $b
671 =item C<gsl_matrix_div_elements($a, $b)> - Divide the elements of matrix $a by the elements of matrix $b
673 =item C<gsl_matrix_scale($a, $x)> - Multiplie the elements of matrix $a by the constant factor $x
675 =item C<gsl_matrix_add_constant($a, $x)> - Add the constant value $x to the elements of the matrix $a
677 =item C<gsl_matrix_add_diagonal($a, $x)> - Add the constant value $x to the elements of the diagonal of the matrix $a
679 =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.
681 =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.
683 =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.
685 =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.
687 =back
689 The following functions are specific to matrices containing complex numbers :
691 =over 1
693 =item C<gsl_matrix_complex_alloc >
695 =item C<gsl_matrix_complex_calloc >
697 =item C<gsl_matrix_complex_alloc_from_block >
699 =item C<gsl_matrix_complex_alloc_from_matrix >
701 =item C<gsl_vector_complex_alloc_row_from_matrix >
703 =item C<gsl_vector_complex_alloc_col_from_matrix >
705 =item C<gsl_matrix_complex_free >
707 =item C<gsl_matrix_complex_submatrix >
709 =item C<gsl_matrix_complex_row >
711 =item C<gsl_matrix_complex_column >
713 =item C<gsl_matrix_complex_diagonal >
715 =item C<gsl_matrix_complex_subdiagonal >
717 =item C<gsl_matrix_complex_superdiagonal >
719 =item C<gsl_matrix_complex_subrow >
721 =item C<gsl_matrix_complex_subcolumn >
723 =item C<gsl_matrix_complex_view_array >
725 =item C<gsl_matrix_complex_view_array_with_tda >
727 =item C<gsl_matrix_complex_view_vector >
729 =item C<gsl_matrix_complex_view_vector_with_tda >
731 =item C<gsl_matrix_complex_const_submatrix >
733 =item C<gsl_matrix_complex_const_row >
735 =item C<gsl_matrix_complex_const_column >
737 =item C<gsl_matrix_complex_const_diagonal >
739 =item C<gsl_matrix_complex_const_subdiagonal >
741 =item C<gsl_matrix_complex_const_superdiagonal >
743 =item C<gsl_matrix_complex_const_subrow >
745 =item C<gsl_matrix_complex_const_subcolumn >
747 =item C<gsl_matrix_complex_const_view_array >
749 =item C<gsl_matrix_complex_const_view_array_with_tda >
751 =item C<gsl_matrix_complex_const_view_vector >
753 =item C<gsl_matrix_complex_const_view_vector_with_tda >
755 =item C<gsl_matrix_complex_get>
757 =item C<gsl_matrix_complex_set>
759 =item C<gsl_matrix_complex_ptr>
761 =item C<gsl_matrix_complex_const_ptr>
763 =item C<gsl_matrix_complex_set_zero >
765 =item C<gsl_matrix_complex_set_identity >
767 =item C<gsl_matrix_complex_set_all >
769 =item C<gsl_matrix_complex_fread >
771 =item C<gsl_matrix_complex_fwrite >
773 =item C<gsl_matrix_complex_fscanf >
775 =item C<gsl_matrix_complex_fprintf >
777 =item C<gsl_matrix_complex_memcpy>
779 =item C<gsl_matrix_complex_swap>
781 =item C<gsl_matrix_complex_swap_rows>
783 =item C<gsl_matrix_complex_swap_columns>
785 =item C<gsl_matrix_complex_swap_rowcol>
787 =item C<gsl_matrix_complex_transpose >
789 =item C<gsl_matrix_complex_transpose_memcpy >
791 =item C<gsl_matrix_complex_isnull >
793 =item C<gsl_matrix_complex_ispos >
795 =item C<gsl_matrix_complex_isneg >
797 =item C<gsl_matrix_complex_add >
799 =item C<gsl_matrix_complex_sub >
801 =item C<gsl_matrix_complex_mul_elements >
803 =item C<gsl_matrix_complex_div_elements >
805 =item C<gsl_matrix_complex_scale >
807 =item C<gsl_matrix_complex_add_constant >
809 =item C<gsl_matrix_complex_add_diagonal >
811 =item C<gsl_matrix_complex_get_row>
813 =item C<gsl_matrix_complex_get_col>
815 =item C<gsl_matrix_complex_set_row>
817 =item C<gsl_matrix_complex_set_col>
819 =back
822 The following functions are the same as the previous enonced ones except that they work with other data types than double.
824 =over 1
826 =item gsl_matrix_const_row gsl_matrix_const_column gsl_matrix_const_diagonal
828 =item gsl_matrix_const_subdiagonal gsl_matrix_const_superdiagonal
830 =item gsl_matrix_const_subrow gsl_matrix_const_subcolumn
832 =item gsl_matrix_const_view_array gsl_matrix_const_view_array_with_tda
834 =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
836 =item gsl_matrix_char_alloc_from_matrix gsl_vector_char_alloc_row_from_matrix gsl_vector_char_alloc_col_from_matrix
838 =item gsl_matrix_char_free gsl_matrix_char_submatrix gsl_matrix_char_row gsl_matrix_char_column
840 =item gsl_matrix_char_diagonal gsl_matrix_char_subdiagonal gsl_matrix_char_superdiagonal
842 =item gsl_matrix_char_subrow gsl_matrix_char_subcolumn gsl_matrix_char_view_array
844 =item gsl_matrix_char_view_array_with_tda gsl_matrix_char_view_vector gsl_matrix_char_view_vector_with_tda
846 =item gsl_matrix_char_const_submatrix gsl_matrix_char_const_row gsl_matrix_char_const_column
848 =item gsl_matrix_char_const_diagonal gsl_matrix_char_const_subdiagonal gsl_matrix_char_const_superdiagonal
850 =item gsl_matrix_char_const_subrow gsl_matrix_char_const_subcolumn gsl_matrix_char_const_view_array
852 =item gsl_matrix_char_const_view_array_with_tda gsl_matrix_char_const_view_vector gsl_matrix_char_const_view_vector_with_tda
854 =item gsl_matrix_char_get gsl_matrix_char_set gsl_matrix_char_ptr gsl_matrix_char_const_ptr
856 =item gsl_matrix_char_set_zero gsl_matrix_char_set_identity
858 =item gsl_matrix_char_set_all gsl_matrix_char_fread
860 =item gsl_matrix_char_fwrite gsl_matrix_char_fscanf gsl_matrix_char_fprintf
862 =item gsl_matrix_char_memcpy gsl_matrix_char_swap
864 =item gsl_matrix_char_swap_rows gsl_matrix_char_swap_columns
866 =item gsl_matrix_char_swap_rowcol gsl_matrix_char_transpose gsl_matrix_char_transpose_memcpy
868 =item gsl_matrix_char_max gsl_matrix_char_min
870 =item gsl_matrix_char_minmax gsl_matrix_char_max_index
872 =item gsl_matrix_char_min_index gsl_matrix_char_minmax_index
874 =item gsl_matrix_char_isnull gsl_matrix_char_ispos gsl_matrix_char_isneg
876 =item gsl_matrix_char_isnonneg gsl_matrix_char_add
878 =item gsl_matrix_char_sub gsl_matrix_char_mul_elements gsl_matrix_char_div_elements
880 =item gsl_matrix_char_scale gsl_matrix_char_add_constant gsl_matrix_char_add_diagonal
882 =item gsl_matrix_int_alloc gsl_matrix_int_calloc gsl_matrix_int_alloc_from_block
884 =item gsl_matrix_int_alloc_from_matrix gsl_vector_int_alloc_row_from_matrix gsl_vector_int_alloc_col_from_matrix
886 =item gsl_matrix_int_free gsl_matrix_int_submatrix gsl_matrix_int_row
888 =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
890 =item gsl_matrix_int_const_submatrix gsl_matrix_int_const_row gsl_matrix_int_const_column
893 =item gsl_matrix_int_const_diagonal gsl_matrix_int_const_subdiagonal gsl_matrix_int_const_superdiagonal
895 =item gsl_matrix_int_const_subrow gsl_matrix_int_const_subcolumn gsl_matrix_int_const_view_array
897 =item gsl_matrix_int_const_view_array_with_tda gsl_matrix_int_const_view_vector gsl_matrix_int_const_view_vector_with_tda
899 =item gsl_matrix_int_get gsl_matrix_int_set
901 =item gsl_matrix_int_ptr gsl_matrix_int_const_ptr
903 =item gsl_matrix_int_set_zero gsl_matrix_int_set_identity gsl_matrix_int_set_all
905 =item gsl_matrix_int_fread gsl_matrix_int_fwrite
907 =item gsl_matrix_int_fscanf gsl_matrix_int_fprintf
909 =item gsl_matrix_int_memcpy gsl_matrix_int_swap
911 =item gsl_matrix_int_swap_rows gsl_matrix_int_swap_columns gsl_matrix_int_swap_rowcol
913 =item gsl_matrix_int_transpose gsl_matrix_int_transpose_memcpy
915 =item gsl_matrix_int_max gsl_matrix_int_min gsl_matrix_int_minmax
917 =item gsl_matrix_int_max_index gsl_matrix_int_min_index
919 =item gsl_matrix_int_minmax_index gsl_matrix_int_isnull
921 =item gsl_matrix_int_ispos gsl_matrix_int_isneg gsl_matrix_int_isnonneg
923 =item gsl_matrix_int_add gsl_matrix_int_sub
925 =item gsl_matrix_int_mul_elements gsl_matrix_int_div_elements gsl_matrix_int_scale
927 =item gsl_matrix_int_add_constant gsl_matrix_int_add_diagonal
929 =back
931 You have to add the functions you want to use inside the qw /put_funtion_here /.
932 You can also write use Math::GSL::Matrix qw/:all/ to use all avaible functions of the module.
933 Other tags are also avaible, here is a complete list of all tags for this module :
935 =over 1
937 =item C<all>
939 =item C<int>
941 =item C<double>
943 =item C<char>
945 =item C<complex>
947 =back
949 For more informations on the functions, we refer you to the GSL offcial documentation: http://www.gnu.org/software/gsl/manual/html_node/
950 Tip : search on google: site:http://www.gnu.org/software/gsl/manual/html_node/name_of_the_function_you_want
952 =head1 EXAMPLES
954 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
956 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.
958 use Math::GSL::Matrix qw/:all/;
959 $m = gsl_matrix_alloc (10,3);
960 for($i = 0; $i< 10; $i++){
961 for($j = 0; $j<3; $j++){
962 gsl_matrix_set($m, $i, $j, 0.23 + 100*$i + $j);
966 for($i = 0; $i< 100; $i++){ # OUT OF RANGE ERROR
967 for($j=0; $j<3; $j++){
968 print "m($i, $j) = " . gsl_matrix_get ($m, $i, $j) . "\n";
971 gsl_matrix_free ($m);
974 use Math::GSL::Matrix qw/:all/;
976 $m = gsl_matrix_alloc (100, 100);
977 $a = gsl_matrix_alloc (100, 100);
979 for($i = 0; $i < 100; $i++){
980 for($j = 0; $j < 100; $j++){
981 gsl_matrix_set ($m, $i, $j, 0.23 + $i + $j);
985 The next program shows how to write a matrix to a file.
987 $out = fopen("test.dat", "wb");
988 gsl_matrix_fwrite ($out, $m);
989 fclose ($out);
991 $in = fopen("test.dat", "rb");
992 gsl_matrix_fread ($in, $a);
993 fclose($in);
995 $k=0;
996 for($i = 0; $i < 100; $i++){
997 for($j = 0; $j < 100; $j++){
998 $mij = gsl_matrix_get ($m, $i, $j);
999 $aij = gsl_matrix_get ($a, $i, $j);
1000 if ($mij != $aij){ $k++ };
1004 gsl_matrix_free($m);
1005 gsl_matrix_free($a);
1007 print "differences = $k (should be zero)\n";
1009 =head1 AUTHOR
1011 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
1013 =head1 COPYRIGHT AND LICENSE
1015 Copyright (C) 2008 Jonathan Leto and Thierry Moisan
1017 This program is free software; you can redistribute it and/or modify it
1018 under the same terms as Perl itself.
1020 =cut