1 %module
"Math::GSL::Linalg"
2 %apply int
*OUTPUT { int
*signum
};
5 #include
"gsl/gsl_linalg.h"
6 #include
"gsl/gsl_permutation.h"
9 %include
"gsl/gsl_linalg.h"
10 %include
"gsl/gsl_permutation.h"
13 @EXPORT_OK
= qw
/$GSL_LINALG_MOD_NONE $GSL_LINALG_MOD_TRANSPOSE $GSL_LINALG_MOD_CONJUGATE
14 gsl_linalg_matmult gsl_linalg_matmult_mod
15 gsl_linalg_exponential_ss
16 gsl_linalg_householder_transform
17 gsl_linalg_complex_householder_transform
18 gsl_linalg_householder_hm
19 gsl_linalg_householder_mh
20 gsl_linalg_householder_hv
21 gsl_linalg_householder_hm1
22 gsl_linalg_complex_householder_hm
23 gsl_linalg_complex_householder_mh
24 gsl_linalg_complex_householder_hv
25 gsl_linalg_hessenberg_decomp
26 gsl_linalg_hessenberg_unpack
27 gsl_linalg_hessenberg_unpack_accum
28 gsl_linalg_hessenberg_set_zero
29 gsl_linalg_hessenberg_submatrix
31 gsl_linalg_hesstri_decomp
33 gsl_linalg_SV_decomp_mod
34 gsl_linalg_SV_decomp_jacobi
44 gsl_linalg_complex_LU_decomp
45 gsl_linalg_complex_LU_solve
46 gsl_linalg_complex_LU_svx
47 gsl_linalg_complex_LU_refine
48 gsl_linalg_complex_LU_invert
49 gsl_linalg_complex_LU_det
50 gsl_linalg_complex_LU_lndet
51 gsl_linalg_complex_LU_sgndet
66 gsl_linalg_QRPT_decomp
67 gsl_linalg_QRPT_decomp2
70 gsl_linalg_QRPT_QRsolve
71 gsl_linalg_QRPT_Rsolve
73 gsl_linalg_QRPT_update
77 gsl_linalg_LQ_lssolve_T
78 gsl_linalg_LQ_Lsolve_T
86 gsl_linalg_PTLQ_decomp
87 gsl_linalg_PTLQ_decomp2
88 gsl_linalg_PTLQ_solve_T
90 gsl_linalg_PTLQ_LQsolve_T
91 gsl_linalg_PTLQ_Lsolve_T
92 gsl_linalg_PTLQ_Lsvx_T
93 gsl_linalg_PTLQ_update
94 gsl_linalg_cholesky_decomp
95 gsl_linalg_cholesky_solve
96 gsl_linalg_cholesky_svx
97 gsl_linalg_cholesky_decomp_unit
98 gsl_linalg_complex_cholesky_decomp
99 gsl_linalg_complex_cholesky_solve
100 gsl_linalg_complex_cholesky_svx
101 gsl_linalg_symmtd_decomp
102 gsl_linalg_symmtd_unpack
103 gsl_linalg_symmtd_unpack_T
104 gsl_linalg_hermtd_decomp
105 gsl_linalg_hermtd_unpack
106 gsl_linalg_hermtd_unpack_T
109 gsl_linalg_solve_symm_tridiag
110 gsl_linalg_solve_tridiag
111 gsl_linalg_solve_symm_cyc_tridiag
112 gsl_linalg_solve_cyc_tridiag
113 gsl_linalg_bidiag_decomp
114 gsl_linalg_bidiag_unpack
115 gsl_linalg_bidiag_unpack2
116 gsl_linalg_bidiag_unpack_B
117 gsl_linalg_balance_matrix
118 gsl_linalg_balance_accum
119 gsl_linalg_balance_columns
121 %EXPORT_TAGS
= ( all
=>[ @EXPORT_OK
] );
127 Math
::GSL
::Linalg
- Functions for solving linear systems
131 use Math
::GSL
::Linalg qw
/:all
/;
136 Here is a list of all the functions included in this module
:
140 =item gsl_linalg_matmult
142 =item gsl_linalg_matmult_mod
144 =item gsl_linalg_exponential_ss
146 =item gsl_linalg_householder_transform
148 =item gsl_linalg_complex_householder_transform
150 =item gsl_linalg_householder_hm
($tau
, $v
, $A
) - This function applies the Householder matrix P defined by the scalar $tau and the vector $v to the left-hand side of the matrix $A. On output the result P A is stored in $A. The function returns
0 if it succeded
, 1 otherwise.
152 =item gsl_linalg_householder_mh
($tau
, $v
, $A
) - This function applies the Householder matrix P defined by the scalar $tau and the vector $v to the right-hand side of the matrix $A. On output the result A P is stored in $A.
154 =item gsl_linalg_householder_hv
($tau
, $v
, $w
) - This function applies the Householder transformation P defined by the scalar $tau and the vector $v to the vector $w. On output the result P w is stored in $w.
156 =item gsl_linalg_householder_hm1
158 =item gsl_linalg_complex_householder_hm
($tau
, $v
, $A
) - Does the same operation than gsl_linalg_householder_hm but with the complex matrix $A
, the complex value $tau and the complex vector $v.
160 =item gsl_linalg_complex_householder_mh
($tau
, $v
, $A
) - Does the same operation than gsl_linalg_householder_mh but with the complex matrix $A
, the complex value $tau and the complex vector $v.
162 =item gsl_linalg_complex_householder_hv
($tau
, $v
, $w
) - Does the same operation than gsl_linalg_householder_hv but with the complex value $tau and the complex vectors $v and $w.
164 =item gsl_linalg_hessenberg_decomp
($A
, $tau
) - This function computes the Hessenberg decomposition of the matrix $A by applying the similarity transformation H
= U^T A U. On output
, H is stored in the upper portion of $A. The information required to construct the matrix U is stored in the lower triangular portion of $A. U is a product of N
- 2 Householder matrices. The Householder vectors are stored in the lower portion of $A
(below the subdiagonal
) and the Householder coefficients are stored in the vector $tau. tau must be of length N. The function returns
0 if it succeded
, 1 otherwise.
166 =item gsl_linalg_hessenberg_unpack
($H
, $tau
, $U
) - This function constructs the orthogonal matrix $U from the information stored in the Hessenberg matrix $H along with the vector $tau. $H and $tau are outputs from gsl_linalg_hessenberg_decomp.
168 =item gsl_linalg_hessenberg_unpack_accum
($H
, $tau
, $V
) - This function is similar to gsl_linalg_hessenberg_unpack
, except it accumulates the matrix U into $V
, so that V'
= VU. The matrix $V must be initialized prior to calling this function. Setting $V to the identity matrix provides the same result as gsl_linalg_hessenberg_unpack. If $H is order N
, then $V must have N columns but may have any number of rows.
170 =item gsl_linalg_hessenberg_set_zero
($H
) - This function sets the lower triangular portion of $H
, below the subdiagonal
, to zero. It is useful for clearing out the Householder vectors after calling gsl_linalg_hessenberg_decomp.
172 =item gsl_linalg_hessenberg_submatrix
174 =item gsl_linalg_hessenberg
176 =item gsl_linalg_hesstri_decomp
($A
, $B
, $U
, $V
, $work
) - This function computes the Hessenberg-Triangular decomposition of the matrix pair
($A
, $B
). On output
, H is stored in $A
, and R is stored in $B. If $U and $V are provided
(they may be null
), the similarity transformations are stored in them. Additional workspace of length N is needed in the vector $work.
178 =item gsl_linalg_SV_decomp
($A
, $V
, $S
, $work
) - This function factorizes the M-by-N matrix $A into the singular value decomposition A
= U S V^T for M
>= N. On output the matrix $A is replaced by U. The diagonal elements of the singular value matrix S are stored in the vector $S. The singular values are non-negative and form a non-increasing sequence from S_1 to S_N. The matrix $V contains the elements of V in untransposed form. To form the product U S V^T it is necessary to take the transpose of V. A workspace of length N is required in vector $work. This routine uses the Golub-Reinsch SVD algorithm.
180 =item gsl_linalg_SV_decomp_mod
($A
, $
X, $V
, $S
, $work
) - This function computes the SVD using the modified Golub-Reinsch algorithm
, which is faster for M
>>N. It requires the vector $work of length N and the N-by-N matrix $
X as additional working space. $A and $V are matrices while $S is a vector.
182 =item gsl_linalg_SV_decomp_jacobi
($A
, $V
, $S
) - This function computes the SVD of the M-by-N matrix $A using one-sided Jacobi orthogonalization for M
>= N. The Jacobi method can compute singular values to higher relative accuracy than Golub-Reinsch algorithms. $V is a matrix while $S is a vector.
184 =item gsl_linalg_SV_solve
($U
, $V
, $S
, $b
, $x
) - This function solves the system A x
= b using the singular value decomposition
($U
, $S
, $V
) of A given by gsl_linalg_SV_decomp. Only non-zero singular values are used in computing the solution. The parts of the solution corresponding to singular values of zero are ignored. Other singular values can be edited out by setting them to zero before calling this function. In the over-determined case where A has more rows than columns the system is solved in the least squares sense
, returning the solution x which minimizes ||A x
- b||_2.
186 =item gsl_linalg_LU_decomp
($a
, $p
) - factorize the matrix $a into the LU decomposition PA
= LU. On output the diagonal and upper triangular part of the input matrix A contain the matrix U. The lower triangular part of the input matrix
(excluding the diagonal
) contains L. The diagonal elements of L are unity
, and are not stored. The function returns two value
, the first is
0 if the operation succeeded
, 1 otherwise
, and the second is the sign of the permutation.
188 =item gsl_linalg_LU_solve
($LU
, $p
, $b
, $x
) - This function solves the square system A x
= b using the LU decomposition of the matrix A into
(LU
, p
) given by gsl_linalg_LU_decomp. $LU is a matrix
, $p a permutation and $b and $x are vectors. The function returns
1 if the operation succeded
, 0 otherwise.
190 =item gsl_linalg_LU_svx
($LU
, $p
, $x
) - This function solves the square system A x
= b in-place using the LU decomposition of A into
(LU
,p
). On input $x should contain the right-hand side b
, which is replaced by the solution on output. $LU is a matrix
, $p a permutation and $x is a vector. The function returns
1 if the operation succeded
, 0 otherwise.
192 =item gsl_linalg_LU_refine
($A
, $LU
, $p
, $b
, $x
, $residual
) - This function apply an iterative improvement to $x
, the solution of $A $x
= $b
, using the LU decomposition of $A into
($LU
,$p
). The initial residual $r
= $A $x
- $b
(where $x and $b are vectors
) is also computed and stored in the vector $residual.
194 =item gsl_linalg_LU_invert
($LU
, $p
, $inverse
) - This function computes the inverse of a matrix A from its LU decomposition stored in the matrix $LU and the permutation $p
, storing the result in the matrix $inverse.
196 =item gsl_linalg_LU_det
($LU
, $signum
) - This function returns the determinant of a matrix A from its LU decomposition stored in the $LU matrix. It needs the integer $signum which is the sign of the permutation returned by gsl_linalg_LU_decomp.
198 =item gsl_linalg_LU_lndet
($LU
) - This function returns the logarithm of the absolute value of the determinant of a matrix A
, from its LU decomposition stored in the $LU matrix.
200 =item gsl_linalg_LU_sgndet
($LU
, $signum
) - This functions computes the sign or phase factor of the determinant of a matrix A
, det
(A
)/|det
(A
)|
, from its LU decomposition
, $LU.
202 =item gsl_linalg_complex_LU_decomp
($A
, $p
) - Does the same operation than gsl_linalg_LU_decomp but on the complex matrix $A.
204 =item gsl_linalg_complex_LU_solve
($LU
, $p
, $b
, $x
) - This functions solve the square system A x
= b using the LU decomposition of A into
($LU
, $p
) given by gsl_linalg_complex_LU_decomp.
206 =item gsl_linalg_complex_LU_svx
($LU
, $p
, $x
) - Does the same operation than gsl_linalg_LU_svx but on the complex matrix $LU and the complex vector $x.
208 =item gsl_linalg_complex_LU_refine
($A
, $LU
, $p
, $b
, $x
, $residual
) - Does the same operation than gsl_linalg_LU_refine but on the complex matrices $A and $LU and with the complex vectors $b
, $x and $residual.
210 =item gsl_linalg_complex_LU_invert
($LU
, $p
, $inverse
) - Does the same operation than gsl_linalg_LU_invert but on the complex matrces $LU and $inverse.
212 =item gsl_linalg_complex_LU_det
($LU
, $signum
) - Does the same operation than gsl_linalg_LU_det but on the complex matrix $LU.
214 =item gsl_linalg_complex_LU_lndet
($LU
) - Does the same operation than gsl_linalg_LU_det but on the complex matrix $LU.
216 =item gsl_linalg_complex_LU_sgndet
($LU
, $signum
) - Does the same operation than gsl_linalg_LU_sgndet but on the complex matrix $LU.
218 =item gsl_linalg_QR_decomp
($a
, $tau
) - factorize the M-by-N matrix A into the QR decomposition A
= Q R. On output the diagonal and upper triangular part of the input matrix $a contain the matrix R. The vector $tau and the columns of the lower triangular part of the matrix $a contain the Householder coefficients and Householder vectors which encode the orthogonal matrix Q. The vector tau must be of length k
= min
(M
,N
).
220 =item gsl_linalg_QR_solve
($QR
, $tau
, $b
, $x
) - This function solves the square system A x
= b using the QR decomposition of A into
(QR
, tau
) given by gsl_linalg_QR_decomp. $QR is matrix
, and $tau
, $b and $x are vectors.
222 =item gsl_linalg_QR_svx
($QR
, $tau
, $x
) - This function solves the square system A x
= b in-place using the QR decomposition of A into the matrix $QR and the vector $tau given by gsl_linalg_QR_decomp. On input
, the vector $x should contain the right-hand side b
, which is replaced by the solution on output.
224 =item gsl_linalg_QR_lssolve
($QR
, $tau
, $b
, $x
, $residual
) - This function finds the least squares solution to the overdetermined system $A $x
= $b where the matrix $A has more rows than columns. The least squares solution minimizes the Euclidean norm of the residual
, ||Ax
- b||.The routine uses the $QR decomposition of $A into
($QR
, $tau
) given by gsl_linalg_QR_decomp. The solution is returned in $x. The residual is computed as a by-product and stored in residual. The function returns
0 if it succeded
, 1 otherwise.
226 =item gsl_linalg_QR_QRsolve
($Q
, $R
, $b
, $x
) - This function solves the system $R $x
= $Q
**T $b for $x. It can be used when the $QR decomposition of a matrix is available in unpacked form as
($Q
, $R
). The function returns
0 if it succeded
, 1 otherwise.
228 =item gsl_linalg_QR_Rsolve
($QR
, $b
, $x
) - This function solves the triangular system R $x
= $b for $x. It may be useful if the product b'
= Q^T b has already been computed using gsl_linalg_QR_QTvec.
230 =item gsl_linalg_QR_Rsvx
($QR
, $x
) - This function solves the triangular system R $x
= b for $x in-place. On input $x should contain the right-hand side b and is replaced by the solution on output. This function may be useful if the product b'
= Q^T b has already been computed using gsl_linalg_QR_QTvec. The function returns
0 if it succeded
, 1 otherwise.
232 =item gsl_linalg_QR_update
($Q
, $R
, $b
, $x
) - This function performs a rank-1 update $w $v
**T of the QR decomposition
($Q
, $R
). The update is given by Q'R'
= Q R
+ w v^T where the output matrices Q' and R' are also orthogonal and right triangular. Note that w is destroyed by the update. The function returns
0 if it succeded
, 1 otherwise.
234 =item gsl_linalg_QR_QTvec
($QR
, $tau
, $v
) - This function applies the matrix Q^T encoded in the decomposition
($QR
,$tau
) to the vector $v
, storing the result Q^T v in $v. The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q^T. The function returns
0 if it succeded
, 1 otherwise.
236 =item gsl_linalg_QR_Qvec
($QR
, $tau
, $v
) - This function applies the matrix Q encoded in the decomposition
($QR
,$tau
) to the vector $v
, storing the result Q v in $v. The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q. The function returns
0 if it succeded
, 1 otherwise.
238 =item gsl_linalg_QR_QTmat
($QR
, $tau
, $A
) - This function applies the matrix Q^T encoded in the decomposition
($QR
,$tau
) to the matrix $A
, storing the result Q^T A in $A. The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q^T. The function returns
0 if it succeded
, 1 otherwise.
240 =item gsl_linalg_QR_unpack
($QR
, $tau
, $Q
, $R
) - This function unpacks the encoded QR decomposition
($QR
,$tau
) into the matrices $Q and $R
, where $Q is M-by-M and $R is M-by-N. The function returns
0 if it succeded
, 1 otherwise.
242 =item gsl_linalg_R_solve
($R
, $b
, $x
) - This function solves the triangular system $R $x
= $b for the N-by-N matrix $R. The function returns
0 if it succeded
, 1 otherwise.
244 =item gsl_linalg_R_svx
($R
, $x
) - This function solves the triangular system $R $x
= b in-place. On input $x should contain the right-hand side b
, which is replaced by the solution on output. The function returns
0 if it succeded
, 1 otherwise.
246 =item gsl_linalg_QRPT_decomp
($A
, $tau
, $p
, $norm
) - This function factorizes the M-by-N matrix $A into the QRP^T decomposition A
= Q R P^T. On output the diagonal and upper triangular part of the input matrix contain the matrix R. The permutation matrix P is stored in the permutation $p. There's two value returned by this function
: the first is
0 if the operation succeeded
, 1 otherwise. The second is sign of the permutation. It has the value
(-1)^n
, where n is the number of interchanges in the permutation. The vector $tau and the columns of the lower triangular part of the matrix $A contain the Householder coefficients and vectors which encode the orthogonal matrix Q. The vector tau must be of length k
=\min
(M
,N
). The matrix Q is related to these components by
, Q
= Q_k ... Q_2 Q_1 where Q_i
= I
- \tau_i v_i v_i^T and v_i is the Householder vector v_i
= (0,...
,1,A
(i
+1,i
),A
(i
+2,i
),...
,A
(m
,i
)). This is the same storage scheme as used by lapack. The vector norm is a workspace of length N used for column pivoting. The algorithm used to perform the decomposition is Householder QR with column pivoting
(Golub
& Van Loan, Matrix Computations, Algorithm 5.4.1).
248 =item gsl_linalg_QRPT_decomp2
($A
, $q
, $r
, $tau
, $p
, $norm
) - This function factorizes the matrix $A into the decomposition A
= Q R P^T without modifying $A itself and storing the output in the separate matrices $q and $r. For the rest
, it operates exactly like gsl_linalg_QRPT_decomp
250 =item gsl_linalg_QRPT_solve
252 =item gsl_linalg_QRPT_svx
254 =item gsl_linalg_QRPT_QRsolve
256 =item gsl_linalg_QRPT_Rsolve
258 =item gsl_linalg_QRPT_Rsvx
260 =item gsl_linalg_QRPT_update
262 =item gsl_linalg_LQ_decomp
264 =item gsl_linalg_LQ_solve_T
266 =item gsl_linalg_LQ_svx_T
268 =item gsl_linalg_LQ_lssolve_T
270 =item gsl_linalg_LQ_Lsolve_T
272 =item gsl_linalg_LQ_Lsvx_T
274 =item gsl_linalg_L_solve_T
276 =item gsl_linalg_LQ_vecQ
278 =item gsl_linalg_LQ_vecQT
280 =item gsl_linalg_LQ_unpack
282 =item gsl_linalg_LQ_update
284 =item gsl_linalg_LQ_LQsolve
286 =item gsl_linalg_PTLQ_decomp
288 =item gsl_linalg_PTLQ_decomp2
290 =item gsl_linalg_PTLQ_solve_T
292 =item gsl_linalg_PTLQ_svx_T
294 =item gsl_linalg_PTLQ_LQsolve_T
296 =item gsl_linalg_PTLQ_Lsolve_T
298 =item gsl_linalg_PTLQ_Lsvx_T
300 =item gsl_linalg_PTLQ_update
302 =item gsl_linalg_cholesky_decomp
($A
) - Factorize the symmetric
, positive-definite square matrix $A into the Cholesky decomposition A
= L L^T and stores it into the matrix $A. The funtcion returns
0 if the operation succeeded
, 0 otherwise.
304 =item gsl_linalg_cholesky_solve
($cholesky
, $b
, $x
) - This function solves the system A x
= b using the Cholesky decomposition of A into the matrix $cholesky given by gsl_linalg_cholesky_decomp. $b and $x are vectors. The funtcion returns
0 if the operation succeeded
, 0 otherwise.
306 =item gsl_linalg_cholesky_svx
($cholesky
, $x
) - This function solves the system A x
= b in-place using the Cholesky decomposition of A into the matrix $cholesky given by gsl_linalg_cholesky_decomp. On input the vector $x should contain the right-hand side b
, which is replaced by the solution on output. The funtcion returns
0 if the operation succeeded
, 0 otherwise.
308 =item gsl_linalg_cholesky_decomp_unit
310 =item gsl_linalg_complex_cholesky_decomp
($A
) - Factorize the symmetric
, positive-definite square matrix $A which contains complex numbers into the Cholesky decomposition A
= L L^T and stores it into the matrix $A. The funtcion returns
0 if the operation succeeded
, 0 otherwise.
312 =item gsl_linalg_complex_cholesky_solve
($cholesky
, $b
, $x
) - This function solves the system A x
= b using the Cholesky decomposition of A into the matrix $cholesky given by gsl_linalg_complex_cholesky_decomp. $b and $x are vectors. The funtcion returns
0 if the operation succeeded
, 0 otherwise.
314 =item gsl_linalg_complex_cholesky_svx
($cholesky
, $x
) - This function solves the system A x
= b in-place using the Cholesky decomposition of A into the matrix $cholesky given by gsl_linalg_complex_cholesky_decomp. On input the vector $x should contain the right-hand side b
, which is replaced by the solution on output. The funtcion returns
0 if the operation succeeded
, 0 otherwise.
316 =item gsl_linalg_symmtd_decomp
($A
, $tau
) - This function factorizes the symmetric square matrix $A into the symmetric tridiagonal decomposition Q T Q^T. On output the diagonal and subdiagonal part of the input matrix $A contain the tridiagonal matrix T. The remaining lower triangular part of the input matrix contains the Householder vectors which
, together with the Householder coefficients $tau
, encode the orthogonal matrix Q. This storage scheme is the same as used by lapack. The upper triangular part of $A is not referenced. $tau is a vector.
318 =item gsl_linalg_symmtd_unpack
($A
, $tau
, $Q
, $diag
, $subdiag
) - This function unpacks the encoded symmetric tridiagonal decomposition
($A
, $tau
) obtained from gsl_linalg_symmtd_decomp into the orthogonal matrix $Q
, the vector of diagonal elements $diag and the vector of subdiagonal elements $subdiag.
320 =item gsl_linalg_symmtd_unpack_T
($A
, $diag
, $subdiag
) - This function unpacks the diagonal and subdiagonal of the encoded symmetric tridiagonal decomposition
($A
, $tau
) obtained from gsl_linalg_symmtd_decomp into the vectors $diag and $subdiag.
322 =item gsl_linalg_hermtd_decomp
($A
, $tau
) - This function factorizes the hermitian matrix $A into the symmetric tridiagonal decomposition U T U^T. On output the real parts of the diagonal and subdiagonal part of the input matrix $A contain the tridiagonal matrix T. The remaining lower triangular part of the input matrix contains the Householder vectors which
, together with the Householder coefficients $tau
, encode the orthogonal matrix Q. This storage scheme is the same as used by lapack. The upper triangular part of $A and imaginary parts of the diagonal are not referenced. $A is a complex matrix and $tau a complex vector.
324 =item gsl_linalg_hermtd_unpack
($A
, $tau
, $U
, $diag
, $subdiag
) - This function unpacks the encoded tridiagonal decomposition
($A
, $tau
) obtained from gsl_linalg_hermtd_decomp into the unitary complex matrix $U
, the real vector of diagonal elements $diag and the real vector of subdiagonal elements $subdiag.
326 =item gsl_linalg_hermtd_unpack_T
($A
, $diag
, $subdiag
) - This function unpacks the diagonal and subdiagonal of the encoded tridiagonal decomposition
(A
, tau
) obtained from the gsl_linalg_hermtd_decomp into the real vectors $diag and $subdiag.
328 =item gsl_linalg_HH_solve
($a
, $b
, $x
) - This function solves the system $A $x
= $b directly using Householder transformations where $A is a matrix
, $b and $x vectors. On output the solution is stored in $x and $b is not modified. $A is destroyed by the Householder transformations.
330 =item gsl_linalg_HH_svx
($A
, $x
) - This function solves the system $A $x
= b in-place using Householder transformations where $A is a matrix
, $b is a vector. On input $x should contain the right-hand side b
, which is replaced by the solution on output. The matrix $A is destroyed by the Householder transformations.
332 =item gsl_linalg_solve_symm_tridiag
334 =item gsl_linalg_solve_tridiag
336 =item gsl_linalg_solve_symm_cyc_tridiag
338 =item gsl_linalg_solve_cyc_tridiag
340 =item gsl_linalg_bidiag_decomp
342 =item gsl_linalg_bidiag_unpack
344 =item gsl_linalg_bidiag_unpack2
346 =item gsl_linalg_bidiag_unpack_B
348 =item gsl_linalg_balance_matrix
350 =item gsl_linalg_balance_accum
352 =item gsl_linalg_balance_columns
355 You have to add the functions you want to use inside the qw
/put_funtion_here
/ with spaces between each function. You can also write use Math
::GSL
::Complex qw
/:all
/ to use all avaible functions of the module.
357 For more informations on the functions
, we refer you to the GSL offcial documentation
: L
<http
://www.gnu.org
/software
/gsl
/manual
/html_node
/>
358 Tip
: search on google
: site
:http
://www.gnu.org
/software
/gsl
/manual
/html_node
/ name_of_the_function_you_want
364 This example shows how to compute the determinant of a matrix with the LU decomposition
:
366 use Math
::GSL
::Matrix qw
/:all
/;
367 use Math
::GSL
::Permutation qw
/:all
/;
368 use Math
::GSL
::Linalg qw
/:all
/;
370 my $Matrix
= gsl_matrix_alloc
(4,4);
371 map
{ gsl_matrix_set
($Matrix
, 0, $_
, $_
+1) } (0.
.3);
373 gsl_matrix_set
($Matrix
,1, 0, 2);
374 gsl_matrix_set
($Matrix
, 1, 1, 3);
375 gsl_matrix_set
($Matrix
, 1, 2, 4);
376 gsl_matrix_set
($Matrix
, 1, 3, 1);
378 gsl_matrix_set
($Matrix
, 2, 0, 3);
379 gsl_matrix_set
($Matrix
, 2, 1, 4);
380 gsl_matrix_set
($Matrix
, 2, 2, 1);
381 gsl_matrix_set
($Matrix
, 2, 3, 2);
383 gsl_matrix_set
($Matrix
, 3, 0, 4);
384 gsl_matrix_set
($Matrix
, 3, 1, 1);
385 gsl_matrix_set
($Matrix
, 3, 2, 2);
386 gsl_matrix_set
($Matrix
, 3, 3, 3);
388 my $permutation
= gsl_permutation_alloc
(4);
389 gsl_permutation_init
($permutation
);
390 my
($result
, $signum
) = gsl_linalg_LU_decomp
($Matrix
, $permutation
);
391 my $det
= gsl_linalg_LU_det
($Matrix
, $signum
);
392 print
"The value of the determinant of the matrix is $det \n";
396 Jonathan Leto
<jonathan@leto.net
> and Thierry Moisan
<thierry.moisan@gmail.com
>
398 =head1 COPYRIGHT
AND LICENSE
400 Copyright
(C
) 2008 Jonathan Leto and Thierry Moisan
402 This program is free software
; you can redistribute it and
/or modify it
403 under the same terms as Perl itself.