Adding documentation and testing skeleton for CBLAS
[Math-GSL.git] / Matrix.i
blob97910d060415c9fdc8a8e15a2405b01d0d33f667
1 %module Matrix
3 %apply int *OUTPUT { size_t *imin, size_t *imax, size_t *jmin, size_t *jmax };
5 %apply double *OUTPUT { double * min_out, double * max_out };
6 %{
7 #include "/usr/local/include/gsl/gsl_matrix.h"
8 #include "/usr/local/include/gsl/gsl_vector_double.h"
9 #include "/usr/local/include/gsl/gsl_matrix_double.h"
10 #include "/usr/local/include/gsl/gsl_matrix_int.h"
11 #include "/usr/local/include/gsl/gsl_matrix_complex_double.h"
12 #include "/usr/local/include/gsl/gsl_matrix_char.h"
15 %include "/usr/local/include/gsl/gsl_matrix.h"
16 %include "/usr/local/include/gsl/gsl_vector_double.h"
17 %include "/usr/local/include/gsl/gsl_matrix_double.h"
18 %include "/usr/local/include/gsl/gsl_matrix_int.h"
19 %include "/usr/local/include/gsl/gsl_matrix_complex_double.h"
20 %include "/usr/local/include/gsl/gsl_matrix_char.h"
23 %perlcode %{
25 no warnings 'redefine';
27 @EXPORT_OK = qw/fopen fclose
28 gsl_matrix_alloc gsl_matrix_calloc gsl_matrix_alloc_from_block
29 gsl_matrix_alloc_from_matrix gsl_vector_alloc_row_from_matrix
30 gsl_vector_alloc_col_from_matrix gsl_matrix_free gsl_matrix_submatrix
31 gsl_matrix_row gsl_matrix_column gsl_matrix_diagonal
32 gsl_matrix_subdiagonal gsl_matrix_superdiagonal gsl_matrix_subrow
33 gsl_matrix_subcolumn gsl_matrix_view_array
34 gsl_matrix_view_array_with_tda gsl_matrix_view_vector
35 gsl_matrix_view_vector_with_tda gsl_matrix_const_submatrix
36 gsl_matrix_const_row gsl_matrix_const_column gsl_matrix_const_diagonal
37 gsl_matrix_const_subdiagonal gsl_matrix_const_superdiagonal
38 gsl_matrix_const_subrow gsl_matrix_const_subcolumn
39 gsl_matrix_const_view_array gsl_matrix_const_view_array_with_tda
40 gsl_matrix_const_view_vector gsl_matrix_const_view_vector_with_tda
41 gsl_matrix_get gsl_matrix_set gsl_matrix_ptr gsl_matrix_const_ptr
42 gsl_matrix_set_zero gsl_matrix_set_identity gsl_matrix_set_all
43 gsl_matrix_fread gsl_matrix_fwrite gsl_matrix_fscanf gsl_matrix_fprintf
44 gsl_matrix_memcpy gsl_matrix_swap gsl_matrix_swap_rows
45 gsl_matrix_swap_columns gsl_matrix_swap_rowcol gsl_matrix_transpose
46 gsl_matrix_transpose_memcpy gsl_matrix_max gsl_matrix_min
47 gsl_matrix_minmax gsl_matrix_max_index gsl_matrix_min_index
48 gsl_matrix_minmax_index gsl_matrix_isnull gsl_matrix_ispos
49 gsl_matrix_isneg gsl_matrix_isnonneg gsl_matrix_add gsl_matrix_sub
50 gsl_matrix_mul_elements gsl_matrix_div_elements gsl_matrix_scale
51 gsl_matrix_add_constant gsl_matrix_add_diagonal
52 gsl_matrix_char_alloc gsl_matrix_char_calloc gsl_matrix_char_alloc_from_block
53 gsl_matrix_char_alloc_from_matrix gsl_vector_char_alloc_row_from_matrix gsl_vector_char_alloc_col_from_matrix
54 gsl_matrix_char_free gsl_matrix_char_submatrix
55 gsl_matrix_char_row gsl_matrix_char_column
56 gsl_matrix_char_diagonal gsl_matrix_char_subdiagonal gsl_matrix_char_superdiagonal
57 gsl_matrix_char_subrow gsl_matrix_char_subcolumn gsl_matrix_char_view_array
58 gsl_matrix_char_view_array_with_tda gsl_matrix_char_view_vector gsl_matrix_char_view_vector_with_tda
59 gsl_matrix_char_const_submatrix gsl_matrix_char_const_row gsl_matrix_char_const_column
60 gsl_matrix_char_const_diagonal gsl_matrix_char_const_subdiagonal gsl_matrix_char_const_superdiagonal
61 gsl_matrix_char_const_subrow gsl_matrix_char_const_subcolumn gsl_matrix_char_const_view_array
62 gsl_matrix_char_const_view_array_with_tda gsl_matrix_char_const_view_vector gsl_matrix_char_const_view_vector_with_tda
63 gsl_matrix_char_get gsl_matrix_char_set gsl_matrix_char_ptr gsl_matrix_char_const_ptr
64 gsl_matrix_char_set_zero gsl_matrix_char_set_identity
65 gsl_matrix_char_set_all gsl_matrix_char_fread
66 gsl_matrix_char_fwrite gsl_matrix_char_fscanf gsl_matrix_char_fprintf
67 gsl_matrix_char_memcpy gsl_matrix_char_swap
68 gsl_matrix_char_swap_rows gsl_matrix_char_swap_columns
69 gsl_matrix_char_swap_rowcol gsl_matrix_char_transpose gsl_matrix_char_transpose_memcpy
70 gsl_matrix_char_max gsl_matrix_char_min
71 gsl_matrix_char_minmax gsl_matrix_char_max_index
72 gsl_matrix_char_min_index gsl_matrix_char_minmax_index
73 gsl_matrix_char_isnull gsl_matrix_char_ispos gsl_matrix_char_isneg
74 gsl_matrix_char_isnonneg gsl_matrix_char_add
75 gsl_matrix_char_sub gsl_matrix_char_mul_elements gsl_matrix_char_div_elements
76 gsl_matrix_char_scale gsl_matrix_char_add_constant gsl_matrix_char_add_diagonal
77 gsl_matrix_int_alloc gsl_matrix_int_calloc gsl_matrix_int_alloc_from_block
78 gsl_matrix_int_alloc_from_matrix gsl_vector_int_alloc_row_from_matrix gsl_vector_int_alloc_col_from_matrix
79 gsl_matrix_int_free gsl_matrix_int_submatrix gsl_matrix_int_row
80 gsl_matrix_int_column gsl_matrix_int_diagonal gsl_matrix_int_subdiagonal
81 gsl_matrix_int_superdiagonal gsl_matrix_int_subrow gsl_matrix_int_subcolumn gsl_matrix_int_view_array
82 gsl_matrix_int_view_array_with_tda gsl_matrix_int_view_vector gsl_matrix_int_view_vector_with_tda
83 gsl_matrix_int_const_submatrix gsl_matrix_int_const_row gsl_matrix_int_const_column
84 gsl_matrix_int_const_diagonal gsl_matrix_int_const_subdiagonal gsl_matrix_int_const_superdiagonal
85 gsl_matrix_int_const_subrow gsl_matrix_int_const_subcolumn gsl_matrix_int_const_view_array
86 gsl_matrix_int_const_view_array_with_tda gsl_matrix_int_const_view_vector gsl_matrix_int_const_view_vector_with_tda
87 gsl_matrix_int_get gsl_matrix_int_set
88 gsl_matrix_int_ptr gsl_matrix_int_const_ptr
89 gsl_matrix_int_set_zero gsl_matrix_int_set_identity gsl_matrix_int_set_all
90 gsl_matrix_int_fread gsl_matrix_int_fwrite
91 gsl_matrix_int_fscanf gsl_matrix_int_fprintf
92 gsl_matrix_int_memcpy gsl_matrix_int_swap
93 gsl_matrix_int_swap_rows gsl_matrix_int_swap_columns gsl_matrix_int_swap_rowcol
94 gsl_matrix_int_transpose gsl_matrix_int_transpose_memcpy
95 gsl_matrix_int_max gsl_matrix_int_min gsl_matrix_int_minmax
96 gsl_matrix_int_max_index gsl_matrix_int_min_index
97 gsl_matrix_int_minmax_index gsl_matrix_int_isnull
98 gsl_matrix_int_ispos gsl_matrix_int_isneg gsl_matrix_int_isnonneg
99 gsl_matrix_int_add gsl_matrix_int_sub
100 gsl_matrix_int_mul_elements gsl_matrix_int_div_elements gsl_matrix_int_scale
101 gsl_matrix_int_add_constant gsl_matrix_int_add_diagonal
102 gsl_matrix_get_row gsl_matrix_get_col gsl_matrix_set_row gsl_matrix_set_col
103 gsl_matrix_complex_alloc gsl_matrix_complex_calloc gsl_matrix_complex_alloc_from_block
104 gsl_matrix_complex_alloc_from_matrix gsl_vector_complex_alloc_row_from_matrix gsl_vector_complex_alloc_col_from_matrix
105 gsl_matrix_complex_free gsl_matrix_complex_submatrix gsl_matrix_complex_row
106 gsl_matrix_complex_column gsl_matrix_complex_diagonal gsl_matrix_complex_subdiagonal
107 gsl_matrix_complex_superdiagonal gsl_matrix_complex_subrow gsl_matrix_complex_subcolumn
108 gsl_matrix_complex_view_array gsl_matrix_complex_view_array_with_tda gsl_matrix_complex_view_vector
109 gsl_matrix_complex_view_vector_with_tda gsl_matrix_complex_const_submatrix gsl_matrix_complex_const_row
110 gsl_matrix_complex_const_column gsl_matrix_complex_const_diagonal gsl_matrix_complex_const_subdiagonal
111 gsl_matrix_complex_const_superdiagonal gsl_matrix_complex_const_subrow gsl_matrix_complex_const_subcolumn
112 gsl_matrix_complex_const_view_array gsl_matrix_complex_const_view_array_with_tda gsl_matrix_complex_const_view_vector
113 gsl_matrix_complex_const_view_vector_with_tda gsl_matrix_complex_get gsl_matrix_complex_set
114 gsl_matrix_complex_ptr gsl_matrix_complex_const_ptr gsl_matrix_complex_set_zero
115 gsl_matrix_complex_set_identity gsl_matrix_complex_set_all gsl_matrix_complex_fread
116 gsl_matrix_complex_fwrite gsl_matrix_complex_fscanf gsl_matrix_complex_fprintf
117 gsl_matrix_complex_memcpy gsl_matrix_complex_swap gsl_matrix_complex_swap_rows
118 gsl_matrix_complex_swap_columns gsl_matrix_complex_swap_rowcol gsl_matrix_complex_transpose
119 gsl_matrix_complex_transpose_memcpy gsl_matrix_complex_isnull gsl_matrix_complex_ispos
120 gsl_matrix_complex_isneg gsl_matrix_complex_add gsl_matrix_complex_sub
121 gsl_matrix_complex_mul_elements gsl_matrix_complex_div_elements gsl_matrix_complex_scale
122 gsl_matrix_complex_add_constant gsl_matrix_complex_add_diagonal gsl_matrix_complex_get_row
123 gsl_matrix_complex_get_col gsl_matrix_complex_set_row gsl_matrix_complex_set_col /;
127 %EXPORT_TAGS = ( all => [ @EXPORT_OK ],
128 char => [ qw/
129 gsl_matrix_char_alloc
130 gsl_matrix_char_calloc
131 gsl_matrix_char_alloc_from_block
132 gsl_matrix_char_alloc_from_matrix
133 gsl_vector_char_alloc_row_from_matrix
134 gsl_vector_char_alloc_col_from_matrix
135 gsl_matrix_char_free
136 gsl_matrix_char_submatrix
137 gsl_matrix_char_row
138 gsl_matrix_char_column
139 gsl_matrix_char_diagonal
140 gsl_matrix_char_subdiagonal
141 gsl_matrix_char_superdiagonal
142 gsl_matrix_char_subrow
143 gsl_matrix_char_subcolumn
144 gsl_matrix_char_view_array
145 gsl_matrix_char_view_array_with_tda
146 gsl_matrix_char_view_vector
147 gsl_matrix_char_view_vector_with_tda
148 gsl_matrix_char_const_submatrix
149 gsl_matrix_char_const_row
150 gsl_matrix_char_const_column
151 gsl_matrix_char_const_diagonal
152 gsl_matrix_char_const_subdiagonal
153 gsl_matrix_char_const_superdiagonal
154 gsl_matrix_char_const_subrow
155 gsl_matrix_char_const_subcolumn
156 gsl_matrix_char_const_view_array
157 gsl_matrix_char_const_view_array_with_tda
158 gsl_matrix_char_const_view_vector
159 gsl_matrix_char_const_view_vector_with_tda
160 gsl_matrix_char_get
161 gsl_matrix_char_set
162 gsl_matrix_char_ptr
163 gsl_matrix_char_const_ptr
164 gsl_matrix_char_set_zero
165 gsl_matrix_char_set_identity
166 gsl_matrix_char_set_all
167 gsl_matrix_char_fread
168 gsl_matrix_char_fwrite
169 gsl_matrix_char_fscanf
170 gsl_matrix_char_fprintf
171 gsl_matrix_char_memcpy
172 gsl_matrix_char_swap
173 gsl_matrix_char_swap_rows
174 gsl_matrix_char_swap_columns
175 gsl_matrix_char_swap_rowcol
176 gsl_matrix_char_transpose
177 gsl_matrix_char_transpose_memcpy
178 gsl_matrix_char_max
179 gsl_matrix_char_min
180 gsl_matrix_char_minmax
181 gsl_matrix_char_max_index
182 gsl_matrix_char_min_index
183 gsl_matrix_char_minmax_index
184 gsl_matrix_char_isnull
185 gsl_matrix_char_ispos
186 gsl_matrix_char_isneg
187 gsl_matrix_char_isnonneg
188 gsl_matrix_char_add
189 gsl_matrix_char_sub
190 gsl_matrix_char_mul_elements
191 gsl_matrix_char_div_elements
192 gsl_matrix_char_scale
193 gsl_matrix_char_add_constant
194 gsl_matrix_char_add_diagonal
197 double => [ qw/
198 gsl_matrix_alloc
199 gsl_matrix_calloc
200 gsl_matrix_alloc_from_block
201 gsl_matrix_alloc_from_matrix
202 gsl_vector_alloc_row_from_matrix
203 gsl_vector_alloc_col_from_matrix
204 gsl_matrix_free
205 gsl_matrix_submatrix
206 gsl_matrix_row
207 gsl_matrix_column
208 gsl_matrix_diagonal
209 gsl_matrix_subdiagonal
210 gsl_matrix_superdiagonal
211 gsl_matrix_subrow
212 gsl_matrix_subcolumn
213 gsl_matrix_view_array
214 gsl_matrix_view_array_with_tda
215 gsl_matrix_view_vector
216 gsl_matrix_view_vector_with_tda
217 gsl_matrix_const_submatrix
218 gsl_matrix_const_row
219 gsl_matrix_const_column
220 gsl_matrix_const_diagonal
221 gsl_matrix_const_subdiagonal
222 gsl_matrix_const_superdiagonal
223 gsl_matrix_const_subrow
224 gsl_matrix_const_subcolumn
225 gsl_matrix_const_view_array
226 gsl_matrix_const_view_array_with_tda
227 gsl_matrix_const_view_vector
228 gsl_matrix_const_view_vector_with_tda
229 gsl_matrix_get
230 gsl_matrix_set
231 gsl_matrix_ptr
232 gsl_matrix_const_ptr
233 gsl_matrix_set_zero
234 gsl_matrix_set_identity
235 gsl_matrix_set_all
236 gsl_matrix_fread
237 gsl_matrix_fwrite
238 gsl_matrix_fscanf
239 gsl_matrix_fprintf
240 gsl_matrix_memcpy
241 gsl_matrix_swap
242 gsl_matrix_swap_rows
243 gsl_matrix_swap_columns
244 gsl_matrix_swap_rowcol
245 gsl_matrix_transpose
246 gsl_matrix_transpose_memcpy
247 gsl_matrix_max
248 gsl_matrix_minmax
249 gsl_matrix_max_index
250 gsl_matrix_min_index
251 gsl_matrix_minmax_index
252 gsl_matrix_isnull
253 gsl_matrix_ispos
254 gsl_matrix_isneg
255 gsl_matrix_isnonneg
256 gsl_matrix_add
257 gsl_matrix_mul_elements
258 gsl_matrix_div_elements
259 gsl_matrix_scale
260 gsl_matrix_add_constant
261 gsl_matrix_add_diagonal
263 int => [ qw/
264 gsl_matrix_int_alloc
265 gsl_matrix_int_alloc_from_matrix
266 gsl_matrix_int_free
267 gsl_matrix_int_column
268 gsl_matrix_int_superdiagonal
269 gsl_matrix_int_view_array_with_tda
270 gsl_matrix_int_const_submatrix
271 gsl_matrix_int_const_diagonal
272 gsl_matrix_int_const_subrow
273 gsl_matrix_int_const_view_array_with_tda
274 gsl_matrix_int_get
275 gsl_matrix_int_ptr
276 gsl_matrix_int_set_zero
277 gsl_matrix_int_fread
278 gsl_matrix_int_fscanf
279 gsl_matrix_int_memcpy
280 gsl_matrix_int_swap_rows
281 gsl_matrix_int_transpose
282 gsl_matrix_int_max
283 gsl_matrix_int_max_index
284 gsl_matrix_int_minmax_index
285 gsl_matrix_int_ispos
286 gsl_matrix_int_add
287 gsl_matrix_int_mul_elements
288 gsl_matrix_int_add_constant
290 complex => [ qw/
291 gsl_matrix_complex_alloc
292 gsl_matrix_complex_calloc
293 gsl_matrix_complex_alloc_from_block
294 gsl_matrix_complex_alloc_from_matrix
295 gsl_vector_complex_alloc_row_from_matrix
296 gsl_vector_complex_alloc_col_from_matrix
297 gsl_matrix_complex_free
298 gsl_matrix_complex_submatrix
299 gsl_matrix_complex_row
300 gsl_matrix_complex_column
301 gsl_matrix_complex_diagonal
302 gsl_matrix_complex_subdiagonal
303 gsl_matrix_complex_superdiagonal
304 gsl_matrix_complex_subrow
305 gsl_matrix_complex_subcolumn
306 gsl_matrix_complex_view_array
307 gsl_matrix_complex_view_array_with_tda
308 gsl_matrix_complex_view_vector
309 gsl_matrix_complex_view_vector_with_tda
310 gsl_matrix_complex_const_submatrix
311 gsl_matrix_complex_const_row
312 gsl_matrix_complex_const_column
313 gsl_matrix_complex_const_diagonal
314 gsl_matrix_complex_const_subdiagonal
315 gsl_matrix_complex_const_superdiagonal
316 gsl_matrix_complex_const_subrow
317 gsl_matrix_complex_const_subcolumn
318 gsl_matrix_complex_const_view_array
319 gsl_matrix_complex_const_view_array_with_tda
320 gsl_matrix_complex_const_view_vector
321 gsl_matrix_complex_const_view_vector_with_tda
322 gsl_matrix_complex_get
323 gsl_matrix_complex_set
324 gsl_matrix_complex_ptr
325 gsl_matrix_complex_const_ptr
326 gsl_matrix_complex_set_zero gsl_matrix_complex_set_identity
327 gsl_matrix_complex_set_all
328 gsl_matrix_complex_fread
329 gsl_matrix_complex_fwrite
330 gsl_matrix_complex_fscanf
331 gsl_matrix_complex_fprintf
332 gsl_matrix_complex_memcpygsl_matrix_complex_swap
333 gsl_matrix_complex_swap_rows
334 gsl_matrix_complex_swap_columns
335 gsl_matrix_complex_swap_rowcol
336 gsl_matrix_complex_transpose
337 gsl_matrix_complex_transpose_memcpy
338 gsl_matrix_complex_isnull
339 gsl_matrix_complex_ispos
340 gsl_matrix_complex_isneg
341 gsl_matrix_complex_add
342 gsl_matrix_complex_sub
343 gsl_matrix_complex_mul_elements
344 gsl_matrix_complex_div_elements
345 gsl_matrix_complex_scale
346 gsl_matrix_complex_add_constant
347 gsl_matrix_complex_add_diagonal
348 gsl_matrix_complex_get_row
349 gsl_matrix_complex_get_col
350 gsl_matrix_complex_set_row
351 gsl_matrix_complex_set_col
352 /]);
354 sub new {
355 my ($class, $rows, $cols) = @_;
356 my $this = {};
357 my $matrix;
358 if ( defined $rows && defined $cols &&
359 $rows > 0 && $cols > 0 &&
360 (int $rows == $rows) && (int $cols == $cols)){
362 $matrix = gsl_matrix_alloc($rows,$cols);
363 } else {
364 die __PACKAGE__.'::new($x,$y) - $x and $y must be positive integers';
366 $this->{_matrix} = $matrix;
367 ($this->{_rows}, $this->{_cols}) = ($rows,$cols);
368 bless $this, $class;
370 sub raw { (shift)->{_matrix} }
371 sub rows { (shift)->{_rows} }
372 sub cols { (shift)->{_cols} }
374 __END__
376 =head1 NAME
378 Math::GSL::Matrix - Mathematical functions concerning Matrices
380 =head1 SYPNOPSIS
382 use Math::GSL::Matrix qw/:all/;
384 =head1 DESCRIPTION
386 Here is a list of all the functions included in this module :
388 =over 1
390 =item C<gsl_matrix_alloc($i, $j)> - Return a gsl_matrix of $i rows and $j columns
392 =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
394 =item C<gsl_matrix_alloc_from_block> -
396 =item C<gsl_matrix_free> -
398 =item C<gsl_matrix_alloc_from_matrix > -
400 =item C<gsl_vector_alloc_row_from_matrix> -
402 =item C<gsl_vector_alloc_col_from_matrix > -
404 =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.
406 =item C<gsl_matrix_row($m , $i)> - Return a vector view of the $i-th row of the matrix $m
408 =item C<gsl_matrix_column($m, $j)> - Return a vector view of the $j-th column of the matrix $m
410 =item C<gsl_matrix_diagonal($m)> - Return a vector view of the diagonal of the vector. The matrix doesn't have to be square.
412 =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.
414 =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.
416 =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.
418 =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.
420 =item C<gsl_matrix_view_array> -
422 =item C<gsl_matrix_view_array_with_tda > -
424 =item C<gsl_matrix_view_vector> -
426 =item C<gsl_matrix_view_vector_with_tda> -
428 =item C<gsl_matrix_const_submatrix> -
430 =item C<gsl_matrix_get($m, $i, $j)> - Return the (i,j)-th element of the matrix $m
432 =item C<gsl_matrix_set($m, $i, $j, $x)> - Set the value of the (i,j)-th element of the matrix $m to $x
434 =item C<gsl_matrix_ptr> -
436 =item C<gsl_matrix_const_ptr> -
438 =item C<gsl_matrix_set_zero($m)> - Set all the elements of the matrix $m to zero
440 =item C<gsl_matrix_set_identity($m)> - Set the elements of the matrix $m to the corresponding elements of the identity matrix
442 =item C<gsl_matrix_set_all($m, $x)> - Set all the elements of the matrix $m to the value $x
444 =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
446 =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
448 =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
450 =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
452 =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.
454 =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.
456 =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.
458 =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.
460 =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.
462 =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.
464 =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.
466 =item C<gsl_matrix_max($m)> - Return the maximum value in the matrix $m
468 =item C<gsl_matrix_min($m)> - Return the minimum value in the matrix $m
470 =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.
472 =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.
474 =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.
476 =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
478 =item C<gsl_matrix_isnull($m)> - Return 1 if all the elements of the matrix $m are zero, 0 otherwise
480 =item C<gsl_matrix_ispos($m)> - Return 1 if all the elements of the matrix $m are strictly positve, 0 otherwise
482 =item C<gsl_matrix_isneg($m)> - Return 1 if all the elements of the matrix $m are strictly negative, 0 otherwise
484 =item C<gsl_matrix_isnonneg($m)> - Return 1 if all the elements of the matrix $m are non-negatuive, 0 otherwise
486 =item C<gsl_matrix_add($a, $b)> - Add the elements of matrix $b to the elements of matrix $a
488 =item C<gsl_matrix_sub($a, $b)> - Subtract the elements of matrix $b from the elements of matrix $a
490 =item C<gsl_matrix_mul_elements($a, $b)> - Multiplie the elements of matrix $a by the elements of matrix $b
492 =item C<gsl_matrix_div_elements($a, $b)> - Divide the elements of matrix $a by the elements of matrix $b
494 =item C<gsl_matrix_scale($a, $x)> - Multiplie the elements of matrix $a by the constant factor $x
496 =item C<gsl_matrix_add_constant($a, $x)> - Add the constant value $x to the elements of the matrix $a
498 =item C<gsl_matrix_add_diagonal($a, $x)> - Add the constant value $x to the elements of the diagonal of the matrix $a
500 =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.
502 =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.
504 =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.
506 =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.
508 =back
510 The following functions are specific to matrices containing complex numbers :
512 =over 1
514 =item C<gsl_matrix_complex_alloc >
516 =item C<gsl_matrix_complex_calloc >
518 =item C<gsl_matrix_complex_alloc_from_block >
520 =item C<gsl_matrix_complex_alloc_from_matrix >
522 =item C<gsl_vector_complex_alloc_row_from_matrix >
524 =item C<gsl_vector_complex_alloc_col_from_matrix >
526 =item C<gsl_matrix_complex_free >
528 =item C<gsl_matrix_complex_submatrix >
530 =item C<gsl_matrix_complex_row >
532 =item C<gsl_matrix_complex_column >
534 =item C<gsl_matrix_complex_diagonal >
536 =item C<gsl_matrix_complex_subdiagonal >
538 =item C<gsl_matrix_complex_superdiagonal >
540 =item C<gsl_matrix_complex_subrow >
542 =item C<gsl_matrix_complex_subcolumn >
544 =item C<gsl_matrix_complex_view_array >
546 =item C<gsl_matrix_complex_view_array_with_tda >
548 =item C<gsl_matrix_complex_view_vector >
550 =item C<gsl_matrix_complex_view_vector_with_tda >
552 =item C<gsl_matrix_complex_const_submatrix >
554 =item C<gsl_matrix_complex_const_row >
556 =item C<gsl_matrix_complex_const_column >
558 =item C<gsl_matrix_complex_const_diagonal >
560 =item C<gsl_matrix_complex_const_subdiagonal >
562 =item C<gsl_matrix_complex_const_superdiagonal >
564 =item C<gsl_matrix_complex_const_subrow >
566 =item C<gsl_matrix_complex_const_subcolumn >
568 =item C<gsl_matrix_complex_const_view_array >
570 =item C<gsl_matrix_complex_const_view_array_with_tda >
572 =item C<gsl_matrix_complex_const_view_vector >
574 =item C<gsl_matrix_complex_const_view_vector_with_tda >
576 =item C<gsl_matrix_complex_get>
578 =item C<gsl_matrix_complex_set>
580 =item C<gsl_matrix_complex_ptr>
582 =item C<gsl_matrix_complex_const_ptr>
584 =item C<gsl_matrix_complex_set_zero >
586 =item C<gsl_matrix_complex_set_identity >
588 =item C<gsl_matrix_complex_set_all >
590 =item C<gsl_matrix_complex_fread >
592 =item C<gsl_matrix_complex_fwrite >
594 =item C<gsl_matrix_complex_fscanf >
596 =item C<gsl_matrix_complex_fprintf >
598 =item C<gsl_matrix_complex_memcpy>
600 =item C<gsl_matrix_complex_swap>
602 =item C<gsl_matrix_complex_swap_rows>
604 =item C<gsl_matrix_complex_swap_columns>
606 =item C<gsl_matrix_complex_swap_rowcol>
608 =item C<gsl_matrix_complex_transpose >
610 =item C<gsl_matrix_complex_transpose_memcpy >
612 =item C<gsl_matrix_complex_isnull >
614 =item C<gsl_matrix_complex_ispos >
616 =item C<gsl_matrix_complex_isneg >
618 =item C<gsl_matrix_complex_add >
620 =item C<gsl_matrix_complex_sub >
622 =item C<gsl_matrix_complex_mul_elements >
624 =item C<gsl_matrix_complex_div_elements >
626 =item C<gsl_matrix_complex_scale >
628 =item C<gsl_matrix_complex_add_constant >
630 =item C<gsl_matrix_complex_add_diagonal >
632 =item C<gsl_matrix_complex_get_row>
634 =item C<gsl_matrix_complex_get_col>
636 =item C<gsl_matrix_complex_set_row>
638 =item C<gsl_matrix_complex_set_col>
640 =back
643 The following functions are the same as the previous enonced ones except that they work with other data types than double.
645 =over 1
647 =item gsl_matrix_const_row gsl_matrix_const_column gsl_matrix_const_diagonal
649 =item gsl_matrix_const_subdiagonal gsl_matrix_const_superdiagonal
651 =item gsl_matrix_const_subrow gsl_matrix_const_subcolumn
653 =item gsl_matrix_const_view_array gsl_matrix_const_view_array_with_tda
655 =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
657 =item gsl_matrix_char_alloc_from_matrix gsl_vector_char_alloc_row_from_matrix gsl_vector_char_alloc_col_from_matrix
659 =item gsl_matrix_char_free gsl_matrix_char_submatrix gsl_matrix_char_row gsl_matrix_char_column
661 =item gsl_matrix_char_diagonal gsl_matrix_char_subdiagonal gsl_matrix_char_superdiagonal
663 =item gsl_matrix_char_subrow gsl_matrix_char_subcolumn gsl_matrix_char_view_array
665 =item gsl_matrix_char_view_array_with_tda gsl_matrix_char_view_vector gsl_matrix_char_view_vector_with_tda
667 =item gsl_matrix_char_const_submatrix gsl_matrix_char_const_row gsl_matrix_char_const_column
669 =item gsl_matrix_char_const_diagonal gsl_matrix_char_const_subdiagonal gsl_matrix_char_const_superdiagonal
671 =item gsl_matrix_char_const_subrow gsl_matrix_char_const_subcolumn gsl_matrix_char_const_view_array
673 =item gsl_matrix_char_const_view_array_with_tda gsl_matrix_char_const_view_vector gsl_matrix_char_const_view_vector_with_tda
675 =item gsl_matrix_char_get gsl_matrix_char_set gsl_matrix_char_ptr gsl_matrix_char_const_ptr
677 =item gsl_matrix_char_set_zero gsl_matrix_char_set_identity
679 =item gsl_matrix_char_set_all gsl_matrix_char_fread
681 =item gsl_matrix_char_fwrite gsl_matrix_char_fscanf gsl_matrix_char_fprintf
683 =item gsl_matrix_char_memcpy gsl_matrix_char_swap
685 =item gsl_matrix_char_swap_rows gsl_matrix_char_swap_columns
687 =item gsl_matrix_char_swap_rowcol gsl_matrix_char_transpose gsl_matrix_char_transpose_memcpy
689 =item gsl_matrix_char_max gsl_matrix_char_min
691 =item gsl_matrix_char_minmax gsl_matrix_char_max_index
693 =item gsl_matrix_char_min_index gsl_matrix_char_minmax_index
695 =item gsl_matrix_char_isnull gsl_matrix_char_ispos gsl_matrix_char_isneg
697 =item gsl_matrix_char_isnonneg gsl_matrix_char_add
699 =item gsl_matrix_char_sub gsl_matrix_char_mul_elements gsl_matrix_char_div_elements
701 =item gsl_matrix_char_scale gsl_matrix_char_add_constant gsl_matrix_char_add_diagonal
703 =item gsl_matrix_int_alloc gsl_matrix_int_calloc gsl_matrix_int_alloc_from_block
705 =item gsl_matrix_int_alloc_from_matrix gsl_vector_int_alloc_row_from_matrix gsl_vector_int_alloc_col_from_matrix
707 =item gsl_matrix_int_free gsl_matrix_int_submatrix gsl_matrix_int_row
709 =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
711 =item gsl_matrix_int_const_submatrix gsl_matrix_int_const_row gsl_matrix_int_const_column
714 =item gsl_matrix_int_const_diagonal gsl_matrix_int_const_subdiagonal gsl_matrix_int_const_superdiagonal
716 =item gsl_matrix_int_const_subrow gsl_matrix_int_const_subcolumn gsl_matrix_int_const_view_array
718 =item gsl_matrix_int_const_view_array_with_tda gsl_matrix_int_const_view_vector gsl_matrix_int_const_view_vector_with_tda
720 =item gsl_matrix_int_get gsl_matrix_int_set
722 =item gsl_matrix_int_ptr gsl_matrix_int_const_ptr
724 =item gsl_matrix_int_set_zero gsl_matrix_int_set_identity gsl_matrix_int_set_all
726 =item gsl_matrix_int_fread gsl_matrix_int_fwrite
728 =item gsl_matrix_int_fscanf gsl_matrix_int_fprintf
730 =item gsl_matrix_int_memcpy gsl_matrix_int_swap
732 =item gsl_matrix_int_swap_rows gsl_matrix_int_swap_columns gsl_matrix_int_swap_rowcol
734 =item gsl_matrix_int_transpose gsl_matrix_int_transpose_memcpy
736 =item gsl_matrix_int_max gsl_matrix_int_min gsl_matrix_int_minmax
738 =item gsl_matrix_int_max_index gsl_matrix_int_min_index
740 =item gsl_matrix_int_minmax_index gsl_matrix_int_isnull
742 =item gsl_matrix_int_ispos gsl_matrix_int_isneg gsl_matrix_int_isnonneg
744 =item gsl_matrix_int_add gsl_matrix_int_sub
746 =item gsl_matrix_int_mul_elements gsl_matrix_int_div_elements gsl_matrix_int_scale
748 =item gsl_matrix_int_add_constant gsl_matrix_int_add_diagonal
750 =back
752 You have to add the functions you want to use inside the qw /put_funtion_here /.
753 You can also write use Math::GSL::PowInt qw/:name_of_tag/ to use all avaible functions of the module.
754 Other tags are also avaible, here is a complete list of all tags for this module :
756 =over 1
758 =item C<all>
760 =item C<int>
762 =item C<double>
764 =item C<char>
766 =item C<complex>
768 =back
770 For more informations on the functions, we refer you to the GSL offcial documentation: http://www.gnu.org/software/gsl/manual/html_node/
771 Tip : search on google: site:http://www.gnu.org/software/gsl/manual/html_node/name_of_the_function_you_want
773 =head1 EXAMPLES
775 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
777 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.
779 use Math::GSL::Matrix qw/:all/;
780 $m = gsl_matrix_alloc (10,3);
781 for($i = 0; $i< 10; $i++){
782 for($j = 0; $j<3; $j++){
783 gsl_matrix_set($m, $i, $j, 0.23 + 100*$i + $j);
787 for($i = 0; $i< 100; $i++){ # OUT OF RANGE ERROR
788 for($j=0; $j<3; $j++){
789 print "m($i, $j) = " . gsl_matrix_get ($m, $i, $j) . "\n";
792 gsl_matrix_free ($m);
795 use Math::GSL::Matrix qw/:all/;
797 $m = gsl_matrix_alloc (100, 100);
798 $a = gsl_matrix_alloc (100, 100);
800 for($i = 0; $i < 100; $i++){
801 for($j = 0; $j < 100; $j++){
802 gsl_matrix_set ($m, $i, $j, 0.23 + $i + $j);
806 The next program shows how to write a matrix to a file.
808 $out = fopen("test.dat", "wb");
809 gsl_matrix_fwrite ($out, $m);
810 fclose ($out);
812 $in = fopen("test.dat", "rb");
813 gsl_matrix_fread ($in, $a);
814 fclose($in);
816 $k=0;
817 for($i = 0; $i < 100; $i++){
818 for($j = 0; $j < 100; $j++){
819 $mij = gsl_matrix_get ($m, $i, $j);
820 $aij = gsl_matrix_get ($a, $i, $j);
821 if ($mij != $aij){ $k++ };
825 gsl_matrix_free($m);
826 gsl_matrix_free($a);
828 print "differences = $k (should be zero)\n";
830 =head1 AUTHOR
832 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
834 =head1 COPYRIGHT AND LICENSE
836 Copyright (C) 2008 Jonathan Leto and Thierry Moisan
838 This program is free software; you can redistribute it and/or modify it
839 under the same terms as Perl itself.
841 =cut