Cleanup Deriv for merge.
[Math-GSL.git] / Eigen.i
blob3cc3cdabf6a002dfff9d94c9ea600c6b28a12d6f
1 %module "Math::GSL::Eigen"
2 %{
3 #include "gsl/gsl_eigen.h"
4 %}
6 %include "gsl/gsl_eigen.h"
8 %perlcode %{
10 @EXPORT_OK = qw/
11 gsl_eigen_symm_alloc gsl_eigen_symm_free
12 gsl_eigen_symm gsl_eigen_symmv_alloc gsl_eigen_symmv_free gsl_eigen_symmv
13 gsl_eigen_herm_alloc gsl_eigen_herm_free gsl_eigen_herm gsl_eigen_hermv_alloc
14 gsl_eigen_hermv_free gsl_eigen_hermv gsl_eigen_francis_alloc gsl_eigen_francis_free
15 gsl_eigen_francis_T gsl_eigen_francis gsl_eigen_francis_Z gsl_eigen_nonsymm_alloc
16 gsl_eigen_nonsymm_free gsl_eigen_nonsymm_params gsl_eigen_nonsymm
17 gsl_eigen_nonsymm_Z gsl_eigen_nonsymmv_alloc gsl_eigen_nonsymmv_free
18 gsl_eigen_nonsymmv gsl_eigen_nonsymmv_Z gsl_eigen_gensymm_alloc
19 gsl_eigen_gensymm_free gsl_eigen_gensymm gsl_eigen_gensymm_standardize
20 gsl_eigen_gensymmv_alloc gsl_eigen_gensymmv_free gsl_eigen_gensymmv
21 gsl_eigen_genherm_alloc gsl_eigen_genherm_free gsl_eigen_genherm
22 gsl_eigen_genherm_standardize gsl_eigen_genhermv_alloc gsl_eigen_genhermv_free
23 gsl_eigen_genhermv gsl_eigen_gen_alloc gsl_eigen_gen_free
24 gsl_eigen_gen_params gsl_eigen_gen gsl_eigen_gen_QZ
25 gsl_eigen_genv_alloc gsl_eigen_genv_free gsl_eigen_genv
26 gsl_eigen_genv_QZ gsl_eigen_symmv_sort gsl_eigen_hermv_sort
27 gsl_eigen_nonsymmv_sort gsl_eigen_gensymmv_sort gsl_eigen_genhermv_sort
28 gsl_eigen_genv_sort gsl_schur_gen_eigvals gsl_schur_solve_equation
29 gsl_schur_solve_equation_z gsl_eigen_jacobi gsl_eigen_invert_jacobi
30 $GSL_EIGEN_SORT_VAL_ASC $GSL_EIGEN_SORT_VAL_DESC
31 $GSL_EIGEN_SORT_ABS_ASC $GSL_EIGEN_SORT_ABS_DESC
33 %EXPORT_TAGS = ( all => [ @EXPORT_OK ] );
35 __END__
37 =head1 NAME
39 Math::GSL::Eigen - Functions for computing eigenvalues and eigenvectors of matrices
41 =head1 SYNOPSIS
43 use Math::GSL::Eigen qw/:all/;
45 =head1 DESCRIPTION
47 Here is a list of all the functions included in this module :
49 =over
51 =item gsl_eigen_symm_alloc($n) - This function returns a workspace for computing eigenvalues of n-by-n real symmetric matrices.
53 =item gsl_eigen_symm_free($w) - This function frees the memory associated with the workspace $w.
55 =item gsl_eigen_symm($A, $eval, $w) - This function computes the eigenvalues of the real symmetric matrix $A. Additional workspace of the appropriate size must be provided in $w. The diagonal and lower triangular part of $A are destroyed during the computation, but the strict upper triangular part is not referenced. The eigenvalues are stored in the vector $eval and are unordered.
57 =item gsl_eigen_symmv_alloc($n) - This function returns a workspace for computing eigenvalues and eigenvectors of n-by-n real symmetric matrices.
59 =item gsl_eigen_symmv_free($w) - This function frees the memory associated with the workspace $w.
61 =item gsl_eigen_symmv($A, $eval, $evec, $w) - This function computes the eigenvalues and eigenvectors of the real symmetric matrix $A. Additional workspace of the appropriate size must be provided in $w. The diagonal and lower triangular part of $A are destroyed during the computation, but the strict upper triangular part is not referenced. The eigenvalues are stored in the vector $eval and are unordered. The corresponding eigenvectors are stored in the columns of the matrix $evec.
63 =item gsl_eigen_herm_alloc($n) - This function returns a workspace for computing eigenvalues of n-by-n complex hermitian matrices.
65 =item gsl_eigen_herm_free($w) - This function frees the memory associated with the workspace $w.
67 =item gsl_eigen_herm($A, $eval, $w) - This function computes the eigenvalues of the complex hermitian matrix $A. Additional workspace of the appropriate size must be provided in $w. The diagonal and lower triangular part of $A are destroyed during the computation, but the strict upper triangular part is not referenced. The imaginary parts of the diagonal are assumed to be zero and are not referenced. The eigenvalues are stored in the vector $eval and are unordered.
69 =item gsl_eigen_hermv_alloc($n) - This function returns a workspace for computing eigenvalues and eigenvectors of n-by-n complex hermitian matrices.
71 =item gsl_eigen_hermv_free($w) - This function frees the memory associated with the workspace $w.
73 =item gsl_eigen_hermv($A, $eval, $evec, $w) - This function computes the eigenvalues and eigenvectors of the complex hermitian matrix $A. Additional workspace of the appropriate size must be provided in $w. The diagonal and lower triangular part of $A are destroyed during the computation, but the strict upper triangular part is not referenced. The imaginary parts of the diagonal are assumed to be zero and are not referenced. The eigenvalues are stored in the vector $eval and are unordered. The corresponding complex eigenvectors are stored in the columns of the matrix $evec.
75 =item gsl_eigen_francis_alloc($n) -
77 =item gsl_eigen_francis_free($w) - This function frees the memory associated with the workspace $w.
79 =item gsl_eigen_francis_T
81 =item gsl_eigen_francis
83 =item gsl_eigen_francis_Z
85 =item gsl_eigen_nonsymm_alloc($n) - This function returns a workspace for computing eigenvalues of n-by-n real nonsymmetric matrices.
87 =item gsl_eigen_nonsymm_free($w) - This function frees the memory associated with the workspace $w.
89 =item gsl_eigen_nonsymm_params($compute_t, $balance, $w) - This function sets some parameters which determine how the eigenvalue problem is solved in subsequent calls to gsl_eigen_nonsymm. If $compute_t is set to 1, the full Schur form T will be computed by gsl_eigen_nonsymm. If it is set to 0, T will not be computed (this is the default setting). If balance is set to 1, a balancing transformation is applied to the matrix prior to computing eigenvalues. This transformation is designed to make the rows and columns of the matrix have comparable norms, and can result in more accurate eigenvalues for matrices whose entries vary widely in magnitude.
91 =item gsl_eigen_nonsymm($A, $eval, $w) - This function computes the eigenvalues of the real nonsymmetric matrix $A and stores them in the vector $eval. If T is desired, it is stored in the upper portion of $A on output. Otherwise, on output, the diagonal of $A will contain the 1-by-1 real eigenvalues and 2-by-2 complex conjugate eigenvalue systems, and the rest of $A is destroyed. In rare cases, this function may fail to find all eigenvalues. If this happens, an error code is returned (1) and the number of converged eigenvalues is stored in $w->{n_evals}. The converged eigenvalues are stored in the beginning of $eval.
93 =item gsl_eigen_nonsymm_Z($A, $eval, $Z, $w) - This function is identical to gsl_eigen_nonsymm except it also computes the Schur vectors and stores them into the $Z matrix.
95 =item gsl_eigen_nonsymmv_alloc($n) - This function allocates a workspace for computing eigenvalues and eigenvectors of n-by-n real nonsymmetric matrices.
97 =item gsl_eigen_nonsymmv_free($w) - This function frees the memory associated with the workspace $w.
99 =item gsl_eigen_nonsymmv($A, $eval, $evec, $w) - This function computes eigenvalues and right eigenvectors of the n-by-n real nonsymmetric matrix $A. It first calls gsl_eigen_nonsymm to compute the eigenvalues, Schur form T, and Schur vectors. Then it finds eigenvectors of T and backtransforms them using the Schur vectors. The Schur vectors are destroyed in the process, but can be saved by using gsl_eigen_nonsymmv_Z. The computed eigenvectors are normalized to have unit magnitude. On output, the upper portion of $A contains the Schur form T. If gsl_eigen_nonsymm fails, no eigenvectors are computed, and an error code is returned (1). $eval is a complex vector and $evec is a complex matrix.
101 =item gsl_eigen_nonsymmv_Z($A, $eval, $evec, $Z, $w) - This function is identical to gsl_eigen_nonsymmv except it also saves the Schur vectors into the $Z matrix.
103 =item gsl_eigen_gensymm_alloc($n) - This function allocates a workspace for computing eigenvalues of n-by-n real generalized symmetric-definite eigensystems.
105 =item gsl_eigen_gensymm_free($w) - This function frees the memory associated with the workspace $w.
107 =item gsl_eigen_gensymm($A, $B, $eval, $w) - This function computes the eigenvalues of the real generalized symmetric-definite matrix pair ($A, $B), and stores them in the vector $eval. On output, $B contains its Cholesky decomposition and $A is destroyed.
109 =item gsl_eigen_gensymm_standardize
111 =item gsl_eigen_gensymmv_alloc($n) - This function allocates a workspace for computing eigenvalues and eigenvectors of n-by-n real generalized symmetric-definite eigensystems.
113 =item gsl_eigen_gensymmv_free($w) - This function frees the memory associated with the workspace $w.
115 =item gsl_eigen_gensymmv($A, $B, $eval, $evec, $w) - This function computes the eigenvalues and eigenvectors of the real generalized symmetric-definite matrix pair ($A, $B), and stores them in $eval vector and $evec matrix respectively. The computed eigenvectors are normalized to have unit magnitude. On output, $B contains its Cholesky decomposition and A is destroyed.
117 =item gsl_eigen_genherm_alloc($n) - This function allocates a workspace for computing eigenvalues of n-by-n complex generalized hermitian-definite eigensystems.
119 =item gsl_eigen_genherm_free($w) - This function frees the memory associated with the workspace $w.
121 =item gsl_eigen_genherm($A, $B, $eval, $w) - This function computes the eigenvalues of the complex generalized hermitian-definite matrix pair ($A, $B), and stores them in the $eval vector. On output, $B contains its Cholesky decomposition and $A is destroyed.
123 =item gsl_eigen_genherm_standardize
125 =item gsl_eigen_genhermv_alloc($n) - This function allocates a workspace for computing eigenvalues and eigenvectors of n-by-n complex generalized hermitian-definite eigensystems.
127 =item gsl_eigen_genhermv_free($w) - This function frees the memory associated with the workspace $w.
129 =item gsl_eigen_genhermv($A, $B, $eval, $evec, $w) - This function computes the eigenvalues and eigenvectors of the complex generalized hermitian-definite matrix pair ($A, $B), and stores them in $eval vector and $evec matrix respectively. The computed eigenvectors are normalized to have unit magnitude. On output, $B contains its Cholesky decomposition and $A is destroyed.
131 =item gsl_eigen_gen_alloc($n) - This function allocates a workspace for computing eigenvalues of n-by-n real generalized nonsymmetric eigensystems.
133 =item gsl_eigen_gen_free($w) - This function frees the memory associated with the workspace $w.
135 =item gsl_eigen_gen_params($compute_s, $compute_t, $balance, $w) - This function sets some parameters which determine how the eigenvalue problem is solved in subsequent calls to gsl_eigen_gen. If $compute_s is set to 1, the full Schur form S will be computed by gsl_eigen_gen. If it is set to 0, S will not be computed (this is the default setting). S is a quasi upper triangular matrix with 1-by-1 and 2-by-2 blocks on its diagonal. 1-by-1 blocks correspond to real eigenvalues, and 2-by-2 blocks correspond to complex eigenvalues. If $compute_t is set to 1, the full Schur form T will be computed by gsl_eigen_gen. If it is set to 0, T will not be computed (this is the default setting). T is an upper triangular matrix with non-negative elements on its diagonal. Any 2-by-2 blocks in S will correspond to a 2-by-2 diagonal block in T. The $balance parameter is currently ignored, since generalized balancing is not yet implemented.
137 =item gsl_eigen_gen($A, $B, $alpha, $beta, $w) - This function computes the eigenvalues of the real generalized nonsymmetric matrix pair ($A, $B), and stores them as pairs in ($alpha, $beta), where $alpha is complex and $beta is real, both are vectors. The elements of $beta are normalized to be non-negative. If S is desired, it is stored in $A on output. If T is desired, it is stored in $B on output. The ordering of eigenvalues in ($alpha, $beta) follows the ordering of the diagonal blocks in the Schur forms S and T. In rare cases, this function may fail to find all eigenvalues. If this occurs, an error code is returned (1).
139 =item gsl_eigen_gen_QZ($A, $B, $alpha, $beta, $Q, $Z, $w) - This function is identical to gsl_eigen_gen except it also computes the left and right Schur vectors and stores them into $Q matrix and $Z matrix respectively.
141 =item gsl_eigen_genv_alloc($n) - This function allocates a workspace for computing eigenvalues and eigenvectors of n-by-n real generalized nonsymmetric eigensystems.
143 =item gsl_eigen_genv_free($w) - This function frees the memory associated with the workspace $w.
145 =item gsl_eigen_genv($A, $B, $alpha, $beta, $evec, $w) - This function computes eigenvalues and right eigenvectors of the n-by-n real generalized nonsymmetric matrix pair ($A, $B). The eigenvalues are stored in ($alpha, $beta) where $alpha is a complex vector and $beta a real vector and the eigenvectors are stored in $evec complex matrix. It first calls gsl_eigen_gen to compute the eigenvalues, Schur forms, and Schur vectors. Then it finds eigenvectors of the Schur forms and backtransforms them using the Schur vectors. The Schur vectors are destroyed in the process, but can be saved by using gsl_eigen_genv_QZ. The computed eigenvectors are normalized to have unit magnitude. On output, ($A, $B) contains the generalized Schur form (S, T). If gsl_eigen_gen fails, no eigenvectors are computed, and an error code is returned (1).
147 =item gsl_eigen_genv_QZ($A, $B, $alpha, $beta, $evec, $Q, $Z, $w) - This function is identical to gsl_eigen_genv except it also computes the left and right Schur vectors and stores them into $Q and $Z matrices respectively.
149 =item gsl_eigen_symmv_sort($eval, $evec, $sort_type) - This function simultaneously sorts the eigenvalues stored in the vector $eval and the corresponding real eigenvectors stored in the columns of the matrix $evec according to the value of the parameter $sort_type which is one of the constant included in this module.
151 =item gsl_eigen_hermv_sort($eval, $evec, $sort_type) - This function simultaneously sorts the eigenvalues stored in the vector $eval and the corresponding real eigenvectors stored in the columns of the matrix $evec according to the value of the parameter $sort_type which is one of the constant included in this module.
153 =item gsl_eigen_nonsymmv_sort($eval, $evec, $sort_type) - This function simultaneously sorts the eigenvalues stored in the vector $eval and the corresponding complex eigenvectors stored in the columns of the complex matrix $evec into ascending or descending order according to the value of the parameter $sort_type. Only $GSL_EIGEN_SORT_ABS_ASC and $GSL_EIGEN_SORT_ABS_DESC are supported due to the eigenvalues being complex.
155 =item gsl_eigen_gensymmv_sort($eval, $evec, $sort_type) - This function simultaneously sorts the eigenvalues stored in the vector $eval and the corresponding real eigenvectors stored in the columns of the matrix $evec according to the value of the parameter $sort_type which is one of the constant included in this module.
157 =item gsl_eigen_genhermv_sort($eval, $evec, $sort_type) - This function simultaneously sorts the eigenvalues stored in the vector $eval and the corresponding real eigenvectors stored in the columns of the matrix $evec according to the value of the parameter $sort_type which is one of the constant included in this module.
159 =item gsl_eigen_genv_sort($eval, $evec, $sort_type) - This function simultaneously sorts the eigenvalues stored in the vector $eval and the corresponding complex eigenvectors stored in the columns of the complex matrix $evec into ascending or descending order according to the value of the parameter $sort_type. Only $GSL_EIGEN_SORT_ABS_ASC and $GSL_EIGEN_SORT_ABS_DESC are supported due to the eigenvalues being complex.
161 =item gsl_schur_gen_eigvals
163 =item gsl_schur_solve_equation
165 =item gsl_schur_solve_equation_z
167 =item gsl_eigen_jacobi
169 =item gsl_eigen_invert_jacobi
171 =back
173 This module also includes these constants :
175 =over
177 =item $GSL_EIGEN_SORT_VAL_ASC - ascending order in numerical value
179 =item $GSL_EIGEN_SORT_VAL_DESC - descending order in numerical value
181 =item $GSL_EIGEN_SORT_ABS_ASC - ascending order in magnitude
183 =item $GSL_EIGEN_SORT_ABS_DESC - descending order in magnitude
185 =back
187 For more informations on the functions, we refer you to the GSL offcial documentation:
188 L<http://www.gnu.org/software/gsl/manual/html_node/>
190 Tip : search on google: site:http://www.gnu.org/software/gsl/manual/html_node/ name_of_the_function_you_want
192 =head1 EXAMPLES
194 This example shows how to use the gsl_eigen_symmv functions to find the eigenvalues and eigenvectors of a matrix.
196 use Math::GSL::Vector qw/:all/;
197 use Math::GSL::Matrix qw/:all/;
198 use Math::GSL::Eigen qw/:all/;
199 my $w = gsl_eigen_symmv_alloc(2);
200 my $m = gsl_matrix_alloc(2,2);
201 gsl_matrix_set($m, 0, 0, 2);
202 gsl_matrix_set($m, 0, 1, 1);
203 gsl_matrix_set($m, 1, 0, 1);
204 gsl_matrix_set($m, 1, 1, 2);
205 my $eval = gsl_vector_alloc(2);
206 my $evec = gsl_matrix_alloc(2,2);
207 gsl_eigen_symmv($m, $eval, $evec, $w);
208 gsl_eigen_gensymmv_sort($eval, $evec, $GSL_EIGEN_SORT_ABS_ASC);
209 print "The first eigenvalue is : " . gsl_vector_get($eval, 0) . "\n";
210 print "The second eigenvalue is : " . gsl_vector_get($eval, 1) . "\n";
211 my $x = gsl_matrix_get($evec, 0, 0);
212 my $y = gsl_matrix_get($evec, 0, 1);
213 print "The first eigenvector is [$x, $y] \n";
214 $x = gsl_matrix_get($evec, 1, 0);
215 $y = gsl_matrix_get($evec, 1, 1);
216 print "The second eigenvector is [$x, $y] \n";
218 =head1 AUTHORS
220 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
222 =head1 COPYRIGHT AND LICENSE
224 Copyright (C) 2008 Jonathan Leto and Thierry Moisan
226 This program is free software; you can redistribute it and/or modify it
227 under the same terms as Perl itself.
229 =cut