rgsl/
linear_algebra.rs

1//
2// A rust binding for the GSL library by Guillaume Gomez (guillaume1.gomez@gmail.com)
3//
4
5/*!
6# Linear Algebra
7
8This chapter describes functions for solving linear systems. The library provides linear algebra operations which operate directly on the
9gsl_vector and gsl_matrix objects. These routines use the standard algorithms from Golub & Van Loan’s Matrix Computations with Level-1 and
10Level-2 BLAS calls for efficiency.
11
12## LU Decomposition
13
14A general square matrix A has an LU decomposition into upper and lower triangular matrices,
15
16P A = L U
17where P is a permutation matrix, L is unit lower triangular matrix and U is upper triangular matrix. For square matrices this decomposition
18can be used to convert the linear system A x = b into a pair of triangular systems (L y = P b, U x = y), which can be solved by forward and
19back-substitution. Note that the LU decomposition is valid for singular matrices.
20
21## QR Decomposition
22
23A general rectangular M-by-N matrix A has a QR decomposition into the product of an orthogonal M-by-M square matrix Q (where Q^T Q = I) and
24an M-by-N right-triangular matrix R,
25
26A = Q R
27This decomposition can be used to convert the linear system A x = b into the triangular system R x = Q^T b, which can be solved by back-substitution.
28Another use of the QR decomposition is to compute an orthonormal basis for a set of vectors. The first N columns of Q form an orthonormal
29basis for the range of A, ran(A), when A has full column rank.
30
31## QR Decomposition with Column Pivoting
32
33The QR decomposition can be extended to the rank deficient case by introducing a column permutation P,
34
35A P = Q R
36The first r columns of Q form an orthonormal basis for the range of A for a matrix with column rank r. This decomposition can also be used
37to convert the linear system A x = b into the triangular system R y = Q^T b, x = P y, which can be solved by back-substitution and permutation.
38We denote the QR decomposition with column pivoting by QRP^T since A = Q R P^T.
39
40## Singular Value Decomposition
41
42A general rectangular M-by-N matrix A has a singular value decomposition (SVD) into the product of an M-by-N orthogonal matrix U, an N-by-N
43diagonal matrix of singular values S and the transpose of an N-by-N orthogonal square matrix V,
44
45A = U S V^T
46
47The singular values \sigma_i = S_{ii} are all non-negative and are generally chosen to form a non-increasing sequence \sigma_1 >= \sigma_2 >=
48... >= \sigma_N >= 0.
49
50The singular value decomposition of a matrix has many practical uses. The condition number of the matrix is given by the ratio of the largest
51singular value to the smallest singular value. The presence of a zero singular value indicates that the matrix is singular. The number of
52non-zero singular values indicates the rank of the matrix. In practice singular value decomposition of a rank-deficient matrix will not produce
53exact zeroes for singular values, due to finite numerical precision. Small singular values should be edited by choosing a suitable tolerance.
54
55For a rank-deficient matrix, the null space of A is given by the columns of V corresponding to the zero singular values. Similarly, the range
56of A is given by columns of U corresponding to the non-zero singular values.
57
58Note that the routines here compute the “thin” version of the SVD with U as M-by-N orthogonal matrix. This allows in-place computation and is
59the most commonly-used form in practice. Mathematically, the “full” SVD is defined with U as an M-by-M orthogonal matrix and S as an M-by-N
60diagonal matrix (with additional rows of zeros).
61
62## Cholesky Decomposition
63
64A symmetric, positive definite square matrix A has a Cholesky decomposition into a product of a lower triangular matrix L and its transpose L^T,
65
66A = L L^T
67
68This is sometimes referred to as taking the square-root of a matrix. The Cholesky decomposition can only be carried out when all the eigenvalues
69of the matrix are positive. This decomposition can be used to convert the linear system A x = b into a pair of triangular systems (L y = b,
70L^T x = y), which can be solved by forward and back-substitution.
71
72## Tridiagonal Decomposition of Real Symmetric Matrices
73
74A symmetric matrix A can be factorized by similarity transformations into the form,
75
76A = Q T Q^T
77
78where Q is an orthogonal matrix and T is a symmetric tridiagonal matrix.
79
80## Tridiagonal Decomposition of Hermitian Matrices
81
82A hermitian matrix A can be factorized by similarity transformations into the form,
83
84A = U T U^T
85
86where U is a unitary matrix and T is a real symmetric tridiagonal matrix.
87
88## Hessenberg Decomposition of Real Matrices
89
90A general real matrix A can be decomposed by orthogonal similarity transformations into the form
91
92A = U H U^T
93
94where U is orthogonal and H is an upper Hessenberg matrix, meaning that it has zeros below the first subdiagonal. The Hessenberg reduction
95is the first step in the Schur decomposition for the nonsymmetric eigenvalue problem, but has applications in other areas as well.
96
97## Hessenberg-Triangular Decomposition of Real Matrices
98
99A general real matrix pair (A, B) can be decomposed by orthogonal similarity transformations into the form
100
101A = U H V^T
102B = U R V^T
103
104where U and V are orthogonal, H is an upper Hessenberg matrix, and R is upper triangular. The Hessenberg-Triangular reduction is the first
105step in the generalized Schur decomposition for the generalized eigenvalue problem.
106
107## Bidiagonalization
108
109A general matrix A can be factorized by similarity transformations into the form,
110
111A = U B V^T
112where U and V are orthogonal matrices and B is a N-by-N bidiagonal matrix with non-zero entries only on the diagonal and superdiagonal. The
113size of U is M-by-N and the size of V is N-by-N.
114
115## Householder Transformations
116
117A Householder transformation is a rank-1 modification of the identity matrix which can be used to zero out selected elements of a vector.
118A Householder matrix P takes the form,
119
120P = I - \tau v v^T
121
122where v is a vector (called the Householder vector) and \tau = 2/(v^T v). The functions described in this section use the rank-1 structure
123of the Householder matrix to create and apply Householder transformations efficiently.
124
125## Tridiagonal Systems
126
127The functions described in this section efficiently solve symmetric, non-symmetric and cyclic tridiagonal systems with minimal storage. Note
128that the current implementations of these functions use a variant of Cholesky decomposition, so the tridiagonal matrix must be positive definite.
129For non-positive definite matrices, the functions return the error code ::Sing.
130
131## Balancing
132
133The process of balancing a matrix applies similarity transformations to make the rows and columns have comparable norms. This is useful, for
134example, to reduce roundoff errors in the solution of eigenvalue problems. Balancing a matrix A consists of replacing A with a similar matrix
135
136A' = D^(-1) A D
137
138where D is a diagonal matrix whose entries are powers of the floating point radix.
139
140##14.16 References and Further Reading
141
142Further information on the algorithms described in this section can be found in the following book,
143
144G. H. Golub, C. F. Van Loan, Matrix Computations (3rd Ed, 1996), Johns Hopkins University Press, ISBN 0-8018-5414-8.
145The LAPACK library is described in the following manual,
146
147LAPACK Users’ Guide (Third Edition, 1999), Published by SIAM, ISBN 0-89871-447-8.
148http://www.netlib.org/lapack
149
150The LAPACK source code can be found at the website above, along with an online copy of the users guide.
151
152The Modified Golub-Reinsch algorithm is described in the following paper,
153
154T.F. Chan, “An Improved Algorithm for Computing the Singular Value Decomposition”, ACM Transactions on Mathematical Software, 8 (1982), pp 72–83.
155The Jacobi algorithm for singular value decomposition is described in the following papers,
156
157J.C. Nash, “A one-sided transformation method for the singular value decomposition and algebraic eigenproblem”, Computer Journal, Volume 18, Number
1581 (1975), p 74–76
159J.C. Nash and S. Shlien “Simple algorithms for the partial singular value decomposition”, Computer Journal, Volume 30 (1987), p 268–275.
160James Demmel, Krešimir Veselić, “Jacobi’s Method is more accurate than QR”, Lapack Working Note 15 (LAWN-15), October 1989. Available from netlib,
161http://www.netlib.org/lapack/ in the lawns or lawnspdf directories.
162!*/
163
164use crate::enums;
165use crate::Value;
166use ffi::FFI;
167
168use types::complex::FFFI;
169
170/// Factorise a general N x N matrix A into,
171///
172///  P A = L U
173///
174/// where P is a permutation matrix, L is unit lower triangular and U is upper triangular.
175///
176/// L is stored in the strict lower triangular part of the input matrix. The diagonal elements of L are unity and are not stored.
177///
178/// U is stored in the diagonal and upper triangular part of the input matrix.
179///
180/// P is stored in the permutation p. Column j of P is column k of the identity matrix, where `k = permutation->data[j]`
181///
182/// signum gives the sign of the permutation, (-1)^n, where n is the  number of interchanges in the permutation.
183///
184/// See Golub & Van Loan, Matrix Computations, Algorithm 3.4.1 (Gauss Elimination with Partial Pivoting).
185#[doc(alias = "gsl_linalg_LU_decomp")]
186pub fn LU_decomp(
187    a: &mut crate::MatrixF64,
188    p: &mut ::Permutation,
189    signum: &mut i32,
190) -> Result<(), Value> {
191    let ret = unsafe { sys::gsl_linalg_LU_decomp(a.unwrap_unique(), p.unwrap_unique(), signum) };
192    result_handler!(ret, ())
193}
194
195/// Factorise a general N x N complex matrix A into,
196///
197///   P A = L U
198///
199/// where P is a permutation matrix, L is unit lower triangular and U is upper triangular.
200///
201/// L is stored in the strict lower triangular part of the input matrix. The diagonal elements of L are unity and are not stored.
202///
203/// U is stored in the diagonal and upper triangular part of the input matrix.
204///
205/// P is stored in the permutation p. Column j of P is column k of the identity matrix, where `k = permutation->data[j]`
206///
207/// signum gives the sign of the permutation, (-1)^n, where n is the number of interchanges in the permutation.
208///
209/// See Golub & Van Loan, Matrix Computations, Algorithm 3.4.1 (Gauss Elimination with Partial Pivoting).
210#[doc(alias = "gsl_linalg_complex_LU_decomp")]
211pub fn complex_LU_decomp(
212    a: &mut crate::MatrixComplexF64,
213    p: &mut ::Permutation,
214    signum: &mut i32,
215) -> Result<(), Value> {
216    let ret =
217        unsafe { sys::gsl_linalg_complex_LU_decomp(a.unwrap_unique(), p.unwrap_unique(), signum) };
218    result_handler!(ret, ())
219}
220
221/// This function solves the square system A x = b using the LU decomposition of A into (LU, p) given by LU_decomp or LU_decomp as input.
222#[doc(alias = "gsl_linalg_LU_solve")]
223pub fn LU_solve(
224    lu: &::MatrixF64,
225    p: &::Permutation,
226    b: &::VectorF64,
227    x: &mut crate::VectorF64,
228) -> Result<(), Value> {
229    let ret = unsafe {
230        sys::gsl_linalg_LU_solve(
231            lu.unwrap_shared(),
232            p.unwrap_shared(),
233            b.unwrap_shared(),
234            x.unwrap_unique(),
235        )
236    };
237    result_handler!(ret, ())
238}
239
240/// This function solves the square system A x = b using the LU decomposition of A into (LU, p) given by LU_decomp or LU_decomp as input.
241#[doc(alias = "gsl_linalg_complex_LU_solve")]
242pub fn complex_LU_solve(
243    lu: &::MatrixComplexF64,
244    p: &::Permutation,
245    b: &::VectorComplexF64,
246    x: &mut crate::VectorComplexF64,
247) -> Result<(), Value> {
248    let ret = unsafe {
249        sys::gsl_linalg_complex_LU_solve(
250            lu.unwrap_shared(),
251            p.unwrap_shared(),
252            b.unwrap_shared(),
253            x.unwrap_unique(),
254        )
255    };
256    result_handler!(ret, ())
257}
258
259/// This function solves the square system A x = b in-place using the precomputed LU decomposition of A into (LU,p). On input x should contain
260/// the right-hand side b, which is replaced by the solution on output.
261#[doc(alias = "gsl_linalg_LU_svx")]
262pub fn LU_svx(lu: &::MatrixF64, p: &::Permutation, x: &mut crate::VectorF64) -> Result<(), Value> {
263    let ret =
264        unsafe { sys::gsl_linalg_LU_svx(lu.unwrap_shared(), p.unwrap_shared(), x.unwrap_unique()) };
265    result_handler!(ret, ())
266}
267
268/// This function solves the square system A x = b in-place using the precomputed LU decomposition of A into (LU,p). On input x should contain
269/// the right-hand side b, which is replaced by the solution on output.
270#[doc(alias = "gsl_linalg_complex_LU_svx")]
271pub fn complex_LU_svx(
272    lu: &::MatrixComplexF64,
273    p: &::Permutation,
274    x: &mut crate::VectorComplexF64,
275) -> Result<(), Value> {
276    let ret = unsafe {
277        sys::gsl_linalg_complex_LU_svx(lu.unwrap_shared(), p.unwrap_shared(), x.unwrap_unique())
278    };
279    result_handler!(ret, ())
280}
281
282/// This function applies an iterative improvement to x, the solution of A x = b, from the precomputed LU decomposition of A into (LU,p). The
283/// initial residual r = A x - b is also computed and stored in residual.
284#[doc(alias = "gsl_linalg_LU_refine")]
285pub fn LU_refine(
286    a: &::MatrixF64,
287    lu: &::MatrixF64,
288    p: &::Permutation,
289    b: &::VectorF64,
290    x: &mut crate::VectorF64,
291    residual: &mut crate::VectorF64,
292) -> Result<(), Value> {
293    let ret = unsafe {
294        sys::gsl_linalg_LU_refine(
295            a.unwrap_shared(),
296            lu.unwrap_shared(),
297            p.unwrap_shared(),
298            b.unwrap_shared(),
299            x.unwrap_unique(),
300            residual.unwrap_unique(),
301        )
302    };
303    result_handler!(ret, ())
304}
305
306/// This function applies an iterative improvement to x, the solution of A x = b, from the precomputed LU decomposition of A into (LU,p). The
307/// initial residual r = A x - b is also computed and stored in residual.
308#[doc(alias = "gsl_linalg_complex_LU_refine")]
309pub fn complex_LU_refine(
310    a: &mut crate::MatrixComplexF64,
311    lu: &::MatrixComplexF64,
312    p: &::Permutation,
313    b: &::VectorComplexF64,
314    x: &mut crate::VectorComplexF64,
315    residual: &mut crate::VectorComplexF64,
316) -> Result<(), Value> {
317    let ret = unsafe {
318        sys::gsl_linalg_complex_LU_refine(
319            a.unwrap_unique(),
320            lu.unwrap_shared(),
321            p.unwrap_shared(),
322            b.unwrap_shared(),
323            x.unwrap_unique(),
324            residual.unwrap_unique(),
325        )
326    };
327    result_handler!(ret, ())
328}
329
330/// This function computes the inverse of a matrix A from its LU decomposition (LU,p), storing the result in the matrix inverse. The inverse
331/// is computed by solving the system A x = b for each column of the identity matrix. It is preferable to avoid direct use of the inverse
332/// whenever possible, as the linear solver functions can obtain the same result more efficiently and reliably (consult any introductory
333/// textbook on numerical linear algebra for details).
334#[doc(alias = "gsl_linalg_LU_invert")]
335pub fn LU_invert(
336    lu: &::MatrixF64,
337    p: &::Permutation,
338    inverse: &mut crate::MatrixF64,
339) -> Result<(), Value> {
340    let ret = unsafe {
341        sys::gsl_linalg_LU_invert(
342            lu.unwrap_shared(),
343            p.unwrap_shared(),
344            inverse.unwrap_unique(),
345        )
346    };
347    result_handler!(ret, ())
348}
349
350/// This function computes the inverse of a matrix A from its LU decomposition (LU,p), storing the result in the matrix inverse. The inverse
351/// is computed by solving the system A x = b for each column of the identity matrix. It is preferable to avoid direct use of the inverse
352/// whenever possible, as the linear solver functions can obtain the same result more efficiently and reliably (consult any introductory
353/// textbook on numerical linear algebra for details).
354#[doc(alias = "gsl_linalg_complex_LU_invert")]
355pub fn complex_LU_invert(
356    lu: &::MatrixComplexF64,
357    p: &::Permutation,
358    inverse: &mut crate::MatrixComplexF64,
359) -> Result<(), Value> {
360    let ret = unsafe {
361        sys::gsl_linalg_complex_LU_invert(
362            lu.unwrap_shared(),
363            p.unwrap_shared(),
364            inverse.unwrap_unique(),
365        )
366    };
367    result_handler!(ret, ())
368}
369
370/// This function computes the determinant of a matrix A from its LU decomposition, LU. The determinant is computed as the product of the
371/// diagonal elements of U and the sign of the row permutation signum.
372#[doc(alias = "gsl_linalg_LU_det")]
373pub fn LU_det(lu: &mut crate::MatrixF64, signum: i32) -> f64 {
374    unsafe { sys::gsl_linalg_LU_det(lu.unwrap_unique(), signum) }
375}
376
377/// This function computes the determinant of a matrix A from its LU decomposition, LU. The determinant is computed as the product of the
378/// diagonal elements of U and the sign of the row permutation signum.
379#[doc(alias = "gsl_linalg_complex_LU_det")]
380pub fn complex_LU_det(lu: &mut crate::MatrixComplexF64, signum: i32) -> crate::ComplexF64 {
381    unsafe { sys::gsl_linalg_complex_LU_det(lu.unwrap_unique(), signum).wrap() }
382}
383
384/// These functions compute the logarithm of the absolute value of the determinant of a matrix A, \ln|\det(A)|, from its LU decomposition,
385/// LU. This function may be useful if the direct computation of the determinant would overflow or underflow.
386#[doc(alias = "gsl_linalg_LU_lndet")]
387pub fn LU_lndet(lu: &mut crate::MatrixF64) -> f64 {
388    unsafe { sys::gsl_linalg_LU_lndet(lu.unwrap_unique()) }
389}
390
391/// This function computes the sign or phase factor of the determinant of a matrix A, \det(A)/|\det(A)|, from its LU decomposition, LU.
392#[doc(alias = "gsl_linalg_complex_LU_lndet")]
393pub fn complex_LU_lndet(lu: &mut crate::MatrixComplexF64) -> f64 {
394    unsafe { sys::gsl_linalg_complex_LU_lndet(lu.unwrap_unique()) }
395}
396
397/// This function computes the sign or phase factor of the determinant of a matrix A, \det(A)/|\det(A)|, from its LU decomposition, LU.
398#[doc(alias = "gsl_linalg_LU_sgndet")]
399pub fn LU_sgndet(lu: &mut crate::MatrixF64, signum: i32) -> i32 {
400    unsafe { sys::gsl_linalg_LU_sgndet(lu.unwrap_unique(), signum) }
401}
402
403/// This function computes the sign or phase factor of the determinant of a matrix A, \det(A)/|\det(A)|, from its LU decomposition, LU.
404#[doc(alias = "gsl_linalg_complex_LU_sgndet")]
405pub fn complex_LU_sgndet(lu: &mut crate::MatrixComplexF64, signum: i32) -> crate::ComplexF64 {
406    unsafe { sys::gsl_linalg_complex_LU_sgndet(lu.unwrap_unique(), signum).wrap() }
407}
408
409/// This function factorizes the M-by-N matrix A into the QR decomposition A = Q R. On output the diagonal and upper triangular part of the
410/// input matrix contain the matrix R. The vector tau and the columns of the lower triangular part of the matrix A contain the Householder
411/// coefficients and Householder vectors which encode the orthogonal matrix Q. The vector tau must be of length k=\min(M,N). The matrix Q
412/// 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 =
413/// (0,...,1,A(i+1,i),A(i+2,i),...,A(m,i)). This is the same storage scheme as used by LAPACK.
414///
415/// The algorithm used to perform the decomposition is Householder QR (Golub & Van Loan, Matrix Computations, Algorithm 5.2.1).
416#[doc(alias = "gsl_linalg_QR_decomp")]
417pub fn QR_decomp(a: &mut crate::MatrixF64, tau: &mut crate::VectorF64) -> Result<(), Value> {
418    let ret = unsafe { sys::gsl_linalg_QR_decomp(a.unwrap_unique(), tau.unwrap_unique()) };
419    result_handler!(ret, ())
420}
421
422/// This function solves the square system A x = b using the QR decomposition of A held in (QR, tau) which must have been computed previously
423/// with gsl_linalg_QR_decomp. The least-squares solution for rectangular systems can be found using QR_lssolve.
424#[doc(alias = "gsl_linalg_QR_solve")]
425pub fn QR_solve(
426    qr: &::MatrixF64,
427    tau: &::VectorF64,
428    b: &::VectorF64,
429    x: &mut crate::VectorF64,
430) -> Result<(), Value> {
431    let ret = unsafe {
432        sys::gsl_linalg_QR_solve(
433            qr.unwrap_shared(),
434            tau.unwrap_shared(),
435            b.unwrap_shared(),
436            x.unwrap_unique(),
437        )
438    };
439    result_handler!(ret, ())
440}
441
442/// This function solves the square system A x = b in-place using the QR decomposition of A held in (QR,tau) which must have been computed
443/// previously by gsl_linalg_QR_decomp. On input x should contain the right-hand side b, which is replaced by the solution on output.
444#[doc(alias = "gsl_linalg_QR_svx")]
445pub fn QR_svx(qr: &::MatrixF64, tau: &::VectorF64, x: &mut crate::VectorF64) -> Result<(), Value> {
446    let ret = unsafe {
447        sys::gsl_linalg_QR_svx(qr.unwrap_shared(), tau.unwrap_shared(), x.unwrap_unique())
448    };
449    result_handler!(ret, ())
450}
451
452/// This function finds the least squares solution to the overdetermined system A x = b where the matrix A has more rows than columns. The
453/// least squares solution minimizes the Euclidean norm of the residual, ||Ax - b||.The routine requires as input the QR decomposition of
454/// 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
455/// residual.
456#[doc(alias = "gsl_linalg_QR_lssolve")]
457pub fn QR_lssolve(
458    qr: &::MatrixF64,
459    tau: &::VectorF64,
460    b: &::VectorF64,
461    x: &mut crate::VectorF64,
462    residual: &mut crate::VectorF64,
463) -> Result<(), Value> {
464    let ret = unsafe {
465        sys::gsl_linalg_QR_lssolve(
466            qr.unwrap_shared(),
467            tau.unwrap_shared(),
468            b.unwrap_shared(),
469            x.unwrap_unique(),
470            residual.unwrap_unique(),
471        )
472    };
473    result_handler!(ret, ())
474}
475
476/// 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
477/// multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q^T.
478#[doc(alias = "gsl_linalg_QR_QTvec")]
479pub fn QR_QTvec(
480    qr: &::MatrixF64,
481    tau: &::VectorF64,
482    v: &mut crate::VectorF64,
483) -> Result<(), Value> {
484    let ret = unsafe {
485        sys::gsl_linalg_QR_QTvec(qr.unwrap_shared(), tau.unwrap_shared(), v.unwrap_unique())
486    };
487    result_handler!(ret, ())
488}
489
490/// 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
491/// multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q.
492#[doc(alias = "gsl_linalg_QR_Qvec")]
493pub fn QR_Qvec(qr: &::MatrixF64, tau: &::VectorF64, v: &mut crate::VectorF64) -> Result<(), Value> {
494    let ret = unsafe {
495        sys::gsl_linalg_QR_Qvec(qr.unwrap_shared(), tau.unwrap_shared(), v.unwrap_unique())
496    };
497    result_handler!(ret, ())
498}
499
500/// 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
501/// multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q^T.
502#[doc(alias = "gsl_linalg_QR_QTmat")]
503pub fn QR_QTmat(
504    qr: &::MatrixF64,
505    tau: &::VectorF64,
506    v: &mut crate::MatrixF64,
507) -> Result<(), Value> {
508    let ret = unsafe {
509        sys::gsl_linalg_QR_QTmat(qr.unwrap_shared(), tau.unwrap_shared(), v.unwrap_unique())
510    };
511    result_handler!(ret, ())
512}
513
514/// 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
515/// gsl_linalg_QR_QTvec.
516#[doc(alias = "gsl_linalg_QR_Rsolve")]
517pub fn QR_Rsolve(qr: &::MatrixF64, b: &::VectorF64, x: &mut crate::VectorF64) -> Result<(), Value> {
518    let ret = unsafe {
519        sys::gsl_linalg_QR_Rsolve(qr.unwrap_shared(), b.unwrap_shared(), x.unwrap_unique())
520    };
521    result_handler!(ret, ())
522}
523
524/// 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
525/// 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.
526#[doc(alias = "gsl_linalg_QR_Rsvx")]
527pub fn QR_Rsvx(qr: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
528    let ret = unsafe { sys::gsl_linalg_QR_Rsvx(qr.unwrap_shared(), x.unwrap_unique()) };
529    result_handler!(ret, ())
530}
531
532/// 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.
533#[doc(alias = "gsl_linalg_QR_unpack")]
534pub fn QR_unpack(
535    qr: &::MatrixF64,
536    tau: &::VectorF64,
537    q: &mut crate::MatrixF64,
538    r: &mut crate::MatrixF64,
539) -> Result<(), Value> {
540    let ret = unsafe {
541        sys::gsl_linalg_QR_unpack(
542            qr.unwrap_shared(),
543            tau.unwrap_shared(),
544            q.unwrap_unique(),
545            r.unwrap_unique(),
546        )
547    };
548    result_handler!(ret, ())
549}
550
551/// 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
552/// form as (Q, R).
553#[doc(alias = "gsl_linalg_QR_QRsolve")]
554pub fn QR_QRsolve(
555    q: &mut crate::MatrixF64,
556    r: &mut crate::MatrixF64,
557    b: &::VectorF64,
558    x: &mut crate::VectorF64,
559) -> Result<(), Value> {
560    let ret = unsafe {
561        sys::gsl_linalg_QR_QRsolve(
562            q.unwrap_unique(),
563            r.unwrap_unique(),
564            b.unwrap_shared(),
565            x.unwrap_unique(),
566        )
567    };
568    result_handler!(ret, ())
569}
570
571/// 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
572/// output matrices Q' and R' are also orthogonal and right triangular. Note that w is destroyed by the update.
573#[doc(alias = "gsl_linalg_QR_update")]
574pub fn QR_update(
575    q: &mut crate::MatrixF64,
576    r: &mut crate::MatrixF64,
577    mut w: crate::VectorF64,
578    v: &::VectorF64,
579) -> Result<(), Value> {
580    let ret = unsafe {
581        sys::gsl_linalg_QR_update(
582            q.unwrap_unique(),
583            r.unwrap_unique(),
584            w.unwrap_unique(),
585            v.unwrap_shared(),
586        )
587    };
588    result_handler!(ret, ())
589}
590
591/// This function solves the triangular system R x = b for the N-by-N matrix R.
592#[doc(alias = "gsl_linalg_R_solve")]
593pub fn R_solve(r: &::MatrixF64, b: &::VectorF64, x: &mut crate::VectorF64) -> Result<(), Value> {
594    let ret =
595        unsafe { sys::gsl_linalg_R_solve(r.unwrap_shared(), b.unwrap_shared(), x.unwrap_unique()) };
596    result_handler!(ret, ())
597}
598
599/// 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
600/// the solution on output.
601#[doc(alias = "gsl_linalg_R_svx")]
602pub fn R_svx(r: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
603    let ret = unsafe { sys::gsl_linalg_R_svx(r.unwrap_shared(), x.unwrap_unique()) };
604    result_handler!(ret, ())
605}
606
607/// 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
608/// of the input matrix contain the matrix R. The permutation matrix P is stored in the permutation p. The sign of the permutation is given
609/// by signum. 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
610/// triangular part of the matrix A contain the Householder coefficients and vectors which encode the orthogonal matrix Q. The vector tau must
611/// 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
612/// 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
613/// a workspace of length N used for column pivoting.
614///
615/// The algorithm used to perform the decomposition is Householder QR with column pivoting (Golub & Van Loan, Matrix Computations, Algorithm 5.4.1).
616#[doc(alias = "gsl_linalg_QRPT_decomp")]
617pub fn QRPT_decomp(
618    a: &mut crate::MatrixF64,
619    tau: &mut crate::VectorF64,
620    p: &mut ::Permutation,
621    signum: &mut i32,
622    norm: &mut crate::VectorF64,
623) -> Result<(), Value> {
624    let ret = unsafe {
625        sys::gsl_linalg_QRPT_decomp(
626            a.unwrap_unique(),
627            tau.unwrap_unique(),
628            p.unwrap_unique(),
629            signum,
630            norm.unwrap_unique(),
631        )
632    };
633    result_handler!(ret, ())
634}
635
636/// 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
637/// matrices q and r.
638#[doc(alias = "gsl_linalg_QRPT_decomp2")]
639pub fn QRPT_decomp2(
640    a: &::MatrixF64,
641    q: &mut crate::MatrixF64,
642    r: &mut crate::MatrixF64,
643    tau: &mut crate::VectorF64,
644    p: &mut ::Permutation,
645    signum: &mut i32,
646    norm: &mut crate::VectorF64,
647) -> Result<(), Value> {
648    let ret = unsafe {
649        sys::gsl_linalg_QRPT_decomp2(
650            a.unwrap_shared(),
651            q.unwrap_unique(),
652            r.unwrap_unique(),
653            tau.unwrap_unique(),
654            p.unwrap_unique(),
655            signum,
656            norm.unwrap_unique(),
657        )
658    };
659    result_handler!(ret, ())
660}
661
662/// This function solves the square system A x = b using the QRP^T decomposition of A held in (QR, tau, p) which must have been computed previously
663/// by QRPT_decomp.
664#[doc(alias = "gsl_linalg_QRPT_solve")]
665pub fn QRPT_solve(
666    qr: &::MatrixF64,
667    tau: &::VectorF64,
668    p: &::Permutation,
669    b: &::VectorF64,
670    x: &mut crate::VectorF64,
671) -> Result<(), Value> {
672    let ret = unsafe {
673        sys::gsl_linalg_QRPT_solve(
674            qr.unwrap_shared(),
675            tau.unwrap_shared(),
676            p.unwrap_shared(),
677            b.unwrap_shared(),
678            x.unwrap_unique(),
679        )
680    };
681    result_handler!(ret, ())
682}
683
684/// This function solves the square system A x = b in-place using the QRP^T decomposition of A held in (QR,tau,p). On input x should contain the
685/// right-hand side b, which is replaced by the solution on output.
686#[doc(alias = "gsl_linalg_QRPT_svx")]
687pub fn QRPT_svx(
688    qr: &::MatrixF64,
689    tau: &::VectorF64,
690    p: &::Permutation,
691    x: &mut crate::VectorF64,
692) -> Result<(), Value> {
693    let ret = unsafe {
694        sys::gsl_linalg_QRPT_svx(
695            qr.unwrap_shared(),
696            tau.unwrap_shared(),
697            p.unwrap_shared(),
698            x.unwrap_unique(),
699        )
700    };
701    result_handler!(ret, ())
702}
703
704/// This function solves the square system R P^T x = Q^T b for x. It can be used when the QR decomposition of a matrix is available in unpacked
705/// form as (Q, R).
706#[doc(alias = "gsl_linalg_QRPT_QRsolve")]
707pub fn QRPT_QRsolve(
708    q: &::MatrixF64,
709    r: &::MatrixF64,
710    p: &::Permutation,
711    b: &::VectorF64,
712    x: &mut crate::VectorF64,
713) -> Result<(), Value> {
714    let ret = unsafe {
715        sys::gsl_linalg_QRPT_QRsolve(
716            q.unwrap_shared(),
717            r.unwrap_shared(),
718            p.unwrap_shared(),
719            b.unwrap_shared(),
720            x.unwrap_unique(),
721        )
722    };
723    result_handler!(ret, ())
724}
725
726/// This function performs a rank-1 update w v^T of the QRP^T decomposition (Q, R, p). The update is given by Q'R' = Q (R + w v^T P) where the
727/// output matrices Q' and R' are also orthogonal and right triangular. Note that w is destroyed by the update. The permutation p is not changed.
728#[doc(alias = "gsl_linalg_QRPT_update")]
729pub fn QRPT_update(
730    q: &mut crate::MatrixF64,
731    r: &mut crate::MatrixF64,
732    p: &::Permutation,
733    w: &mut crate::VectorF64,
734    v: &::VectorF64,
735) -> Result<(), Value> {
736    let ret = unsafe {
737        sys::gsl_linalg_QRPT_update(
738            q.unwrap_unique(),
739            r.unwrap_unique(),
740            p.unwrap_shared(),
741            w.unwrap_unique(),
742            v.unwrap_shared(),
743        )
744    };
745    result_handler!(ret, ())
746}
747
748/// This function solves the triangular system R P^T x = b for the N-by-N matrix R contained in QR.
749#[doc(alias = "gsl_linalg_QRPT_Rsolve")]
750pub fn QRPT_Rsolve(
751    qr: &::MatrixF64,
752    p: &::Permutation,
753    b: &::VectorF64,
754    x: &mut crate::VectorF64,
755) -> Result<(), Value> {
756    let ret = unsafe {
757        sys::gsl_linalg_QRPT_Rsolve(
758            qr.unwrap_shared(),
759            p.unwrap_shared(),
760            b.unwrap_shared(),
761            x.unwrap_unique(),
762        )
763    };
764    result_handler!(ret, ())
765}
766
767/// This function solves the triangular system R P^T x = b in-place for the N-by-N matrix R contained in QR. On input x should contain the
768/// right-hand side b, which is replaced by the solution on output.
769#[doc(alias = "gsl_linalg_QRPT_Rsvx")]
770pub fn QRPT_Rsvx(
771    qr: &::MatrixF64,
772    p: &::Permutation,
773    x: &mut crate::VectorF64,
774) -> Result<(), Value> {
775    let ret = unsafe {
776        sys::gsl_linalg_QRPT_Rsvx(qr.unwrap_shared(), p.unwrap_shared(), x.unwrap_unique())
777    };
778    result_handler!(ret, ())
779}
780
781/// 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
782/// 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
783/// 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
784/// necessary to take the transpose of V. A workspace of length N is required in work.
785///
786/// This routine uses the Golub-Reinsch SVD algorithm.
787#[doc(alias = "gsl_linalg_SV_decomp")]
788pub fn SV_decomp(
789    a: &mut crate::MatrixF64,
790    v: &mut crate::MatrixF64,
791    s: &mut crate::VectorF64,
792    work: &mut crate::VectorF64,
793) -> Result<(), Value> {
794    let ret = unsafe {
795        sys::gsl_linalg_SV_decomp(
796            a.unwrap_unique(),
797            v.unwrap_unique(),
798            s.unwrap_unique(),
799            work.unwrap_unique(),
800        )
801    };
802    result_handler!(ret, ())
803}
804
805/// This function computes the SVD using the modified Golub-Reinsch algorithm, which is faster for M>>N. It requires the vector work of length
806/// N and the N-by-N matrix X as additional working space.
807#[doc(alias = "gsl_linalg_SV_decomp_mod")]
808pub fn SV_decomp_mod(
809    a: &mut crate::MatrixF64,
810    x: &mut crate::MatrixF64,
811    v: &mut crate::MatrixF64,
812    s: &mut crate::VectorF64,
813    work: &mut crate::VectorF64,
814) -> Result<(), Value> {
815    let ret = unsafe {
816        sys::gsl_linalg_SV_decomp_mod(
817            a.unwrap_unique(),
818            x.unwrap_unique(),
819            v.unwrap_unique(),
820            s.unwrap_unique(),
821            work.unwrap_unique(),
822        )
823    };
824    result_handler!(ret, ())
825}
826
827/// 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
828/// singular values to higher relative accuracy than Golub-Reinsch algorithms (see references for details).
829#[doc(alias = "gsl_linalg_SV_decomp_jacobi")]
830pub fn SV_decomp_jacobi(
831    a: &mut crate::MatrixF64,
832    v: &mut crate::MatrixF64,
833    s: &mut crate::VectorF64,
834) -> Result<(), Value> {
835    let ret = unsafe {
836        sys::gsl_linalg_SV_decomp_jacobi(a.unwrap_unique(), v.unwrap_unique(), s.unwrap_unique())
837    };
838    result_handler!(ret, ())
839}
840
841/// This function solves the system A x = b using the singular value decomposition (U, S, V) of A which must have been computed previously
842/// with gsl_linalg_SV_decomp.
843///
844/// Only non-zero singular values are used in computing the solution. The parts of the solution corresponding to singular values of zero are
845/// ignored. Other singular values can be edited out by setting them to zero before calling this function.
846///
847/// In the over-determined case where A has more rows than columns the system is solved in the least squares sense, returning the solution
848/// x which minimizes ||A x - b||_2.
849#[doc(alias = "gsl_linalg_SV_solve")]
850pub fn SV_solve(
851    u: &::MatrixF64,
852    v: &::MatrixF64,
853    s: &::VectorF64,
854    b: &::VectorF64,
855    x: &mut crate::VectorF64,
856) -> Result<(), Value> {
857    let ret = unsafe {
858        sys::gsl_linalg_SV_solve(
859            u.unwrap_shared(),
860            v.unwrap_shared(),
861            s.unwrap_shared(),
862            b.unwrap_shared(),
863            x.unwrap_unique(),
864        )
865    };
866    result_handler!(ret, ())
867}
868
869/// This function computes the statistical leverage values h_i of a matrix A using its singular value decomposition (U, S, V) previously computed
870/// with gsl_linalg_SV_decomp. h_i are the diagonal values of the matrix A (A^T A)^{-1} A^T and depend only on the matrix U which is the input to
871/// this function.
872#[doc(alias = "gsl_linalg_SV_leverage")]
873pub fn SV_leverage(u: &::MatrixF64, h: &mut crate::VectorF64) -> Result<(), Value> {
874    let ret = unsafe { sys::gsl_linalg_SV_leverage(u.unwrap_shared(), h.unwrap_unique()) };
875    result_handler!(ret, ())
876}
877
878/// This function factorizes the symmetric, positive-definite square matrix A into the Cholesky decomposition A = L L^T (or A = L L^H for
879/// the complex case). On input, the values from the diagonal and lower-triangular part of the matrix A are used (the upper triangular part
880/// is ignored). On output the diagonal and lower triangular part of the input matrix A contain the matrix L, while the upper triangular part
881/// of the input matrix is overwritten with L^T (the diagonal terms being identical for both L and L^T). If the matrix is not positive-definite
882/// then the decomposition will fail, returning the error code ::Dom.
883///
884/// When testing whether a matrix is positive-definite, disable the error handler first to avoid triggering an error.
885#[doc(alias = "gsl_linalg_cholesky_decomp")]
886pub fn cholesky_decomp(a: &mut crate::MatrixF64) -> Result<(), Value> {
887    let ret = unsafe { sys::gsl_linalg_cholesky_decomp(a.unwrap_unique()) };
888    result_handler!(ret, ())
889}
890
891/// This function factorizes the symmetric, positive-definite square matrix A into the Cholesky decomposition A = L L^T (or A = L L^H for
892/// the complex case). On input, the values from the diagonal and lower-triangular part of the matrix A are used (the upper triangular part
893/// is ignored). On output the diagonal and lower triangular part of the input matrix A contain the matrix L, while the upper triangular part
894/// of the input matrix is overwritten with L^T (the diagonal terms being identical for both L and L^T). If the matrix is not positive-definite
895/// then the decomposition will fail, returning the error code ::Dom.
896///
897/// When testing whether a matrix is positive-definite, disable the error handler first to avoid triggering an error.
898#[doc(alias = "gsl_linalg_complex_cholesky_decomp")]
899pub fn complex_cholesky_decomp(a: &mut crate::MatrixComplexF64) -> Result<(), Value> {
900    let ret = unsafe { sys::gsl_linalg_complex_cholesky_decomp(a.unwrap_unique()) };
901    result_handler!(ret, ())
902}
903
904/// This function solves the system A x = b using the Cholesky decomposition of A held in the matrix cholesky which must have been previously
905/// computed by gsl_linalg_cholesky_decomp or gsl_linalg_complex_cholesky_decomp.
906#[doc(alias = "gsl_linalg_cholesky_solve")]
907pub fn cholesky_solve(
908    cholesky: &::MatrixF64,
909    b: &::VectorF64,
910    x: &mut crate::VectorF64,
911) -> Result<(), Value> {
912    let ret = unsafe {
913        sys::gsl_linalg_cholesky_solve(
914            cholesky.unwrap_shared(),
915            b.unwrap_shared(),
916            x.unwrap_unique(),
917        )
918    };
919    result_handler!(ret, ())
920}
921
922/// This function solves the system A x = b using the Cholesky decomposition of A held in the matrix cholesky which must have been previously
923/// computed by gsl_linalg_cholesky_decomp or gsl_linalg_complex_cholesky_decomp.
924#[doc(alias = "gsl_linalg_complex_cholesky_solve")]
925pub fn complex_cholesky_solve(
926    cholesky: &::MatrixComplexF64,
927    b: &::VectorComplexF64,
928    x: &mut crate::VectorComplexF64,
929) -> Result<(), Value> {
930    let ret = unsafe {
931        sys::gsl_linalg_complex_cholesky_solve(
932            cholesky.unwrap_shared(),
933            b.unwrap_shared(),
934            x.unwrap_unique(),
935        )
936    };
937    result_handler!(ret, ())
938}
939
940/// This function solves the system A x = b in-place using the Cholesky decomposition of A held in the matrix cholesky which must have been
941/// previously computed by gsl_linalg_cholesky_decomp or gsl_linalg_complex_cholesky_decomp. On input x should contain the right-hand side
942/// b, which is replaced by the solution on output.
943#[doc(alias = "gsl_linalg_cholesky_svx")]
944pub fn cholesky_svx(cholesky: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
945    let ret = unsafe { sys::gsl_linalg_cholesky_svx(cholesky.unwrap_shared(), x.unwrap_unique()) };
946    result_handler!(ret, ())
947}
948
949/// This function solves the system A x = b in-place using the Cholesky decomposition of A held in the matrix cholesky which must have been
950/// previously computed by gsl_linalg_cholesky_decomp or gsl_linalg_complex_cholesky_decomp. On input x should contain the right-hand side
951/// b, which is replaced by the solution on output.
952#[doc(alias = "gsl_linalg_complex_cholesky_svx")]
953pub fn complex_cholesky_svx(
954    cholesky: &::MatrixComplexF64,
955    x: &mut crate::VectorComplexF64,
956) -> Result<(), Value> {
957    let ret = unsafe {
958        sys::gsl_linalg_complex_cholesky_svx(cholesky.unwrap_shared(), x.unwrap_unique())
959    };
960    result_handler!(ret, ())
961}
962
963/// This function computes the inverse of a matrix from its Cholesky decomposition cholesky, which must have been previously computed by
964/// gsl_linalg_cholesky_decomp or gsl_linalg_complex_cholesky_decomp. On output, the inverse is stored in-place in cholesky.
965#[doc(alias = "gsl_linalg_cholesky_invert")]
966pub fn cholesky_invert(cholesky: &mut crate::MatrixF64) -> Result<(), Value> {
967    let ret = unsafe { sys::gsl_linalg_cholesky_invert(cholesky.unwrap_unique()) };
968    result_handler!(ret, ())
969}
970
971/// This function computes the inverse of a matrix from its Cholesky decomposition cholesky, which must have been previously computed by
972/// gsl_linalg_cholesky_decomp or gsl_linalg_complex_cholesky_decomp. On output, the inverse is stored in-place in cholesky.
973#[doc(alias = "gsl_linalg_complex_cholesky_invert")]
974pub fn complex_cholesky_invert(cholesky: &mut crate::MatrixComplexF64) -> Result<(), Value> {
975    let ret = unsafe { sys::gsl_linalg_complex_cholesky_invert(cholesky.unwrap_unique()) };
976    result_handler!(ret, ())
977}
978
979/// This function factorizes the symmetric square matrix A into the symmetric tridiagonal decomposition Q T Q^T. On output the diagonal and
980/// subdiagonal part of the input matrix A contain the tridiagonal matrix T. The remaining lower triangular part of the input matrix contains
981/// the Householder vectors which, together with the Householder coefficients tau, encode the orthogonal matrix Q. This storage scheme is
982/// the same as used by LAPACK. The upper triangular part of A is not referenced.
983#[doc(alias = "gsl_linalg_symmtd_decomp")]
984pub fn symmtd_decomp(a: &mut crate::MatrixF64, tau: &mut crate::VectorF64) -> Result<(), Value> {
985    let ret = unsafe { sys::gsl_linalg_symmtd_decomp(a.unwrap_unique(), tau.unwrap_unique()) };
986    result_handler!(ret, ())
987}
988
989/// This function unpacks the encoded symmetric tridiagonal decomposition (A, tau) obtained from gsl_linalg_symmtd_decomp into the orthogonal
990/// matrix Q, the vector of diagonal elements diag and the vector of subdiagonal elements subdiag.
991#[doc(alias = "gsl_linalg_symmtd_unpack")]
992pub fn symmtd_unpack(
993    a: &::MatrixF64,
994    tau: &::VectorF64,
995    q: &mut crate::MatrixF64,
996    diag: &mut crate::VectorF64,
997    subdiag: &mut crate::VectorF64,
998) -> Result<(), Value> {
999    let ret = unsafe {
1000        sys::gsl_linalg_symmtd_unpack(
1001            a.unwrap_shared(),
1002            tau.unwrap_shared(),
1003            q.unwrap_unique(),
1004            diag.unwrap_unique(),
1005            subdiag.unwrap_unique(),
1006        )
1007    };
1008    result_handler!(ret, ())
1009}
1010
1011/// This function unpacks the diagonal and subdiagonal of the encoded symmetric tridiagonal decomposition (A, tau) obtained from
1012/// gsl_linalg_symmtd_decomp into the vectors diag and subdiag.
1013#[doc(alias = "gsl_linalg_symmtd_unpack_T")]
1014pub fn symmtd_unpack_T(
1015    a: &::MatrixF64,
1016    diag: &mut crate::VectorF64,
1017    subdiag: &mut crate::VectorF64,
1018) -> Result<(), Value> {
1019    let ret = unsafe {
1020        sys::gsl_linalg_symmtd_unpack_T(
1021            a.unwrap_shared(),
1022            diag.unwrap_unique(),
1023            subdiag.unwrap_unique(),
1024        )
1025    };
1026    result_handler!(ret, ())
1027}
1028
1029/// This function factorizes the hermitian matrix A into the symmetric tridiagonal decomposition U T U^T. On output the real parts of the
1030/// diagonal and subdiagonal part of the input matrix A contain the tridiagonal matrix T. The remaining lower triangular part of the input
1031/// matrix contains the Householder vectors which, together with the Householder coefficients tau, encode the unitary matrix U. This storage
1032/// scheme is the same as used by LAPACK. The upper triangular part of A and imaginary parts of the diagonal are not referenced.
1033#[doc(alias = "gsl_linalg_hermtd_decomp")]
1034pub fn hermtd_decomp(
1035    a: &mut crate::MatrixComplexF64,
1036    tau: &mut crate::VectorComplexF64,
1037) -> Result<(), Value> {
1038    let ret = unsafe { sys::gsl_linalg_hermtd_decomp(a.unwrap_unique(), tau.unwrap_unique()) };
1039    result_handler!(ret, ())
1040}
1041
1042/// This function unpacks the encoded tridiagonal decomposition (A, tau) obtained from gsl_linalg_hermtd_decomp into the unitary matrix U,
1043/// the real vector of diagonal elements diag and the real vector of subdiagonal elements subdiag.
1044#[doc(alias = "gsl_linalg_hermtd_unpack")]
1045pub fn hermtd_unpack(
1046    a: &::MatrixComplexF64,
1047    tau: &::VectorComplexF64,
1048    u: &mut crate::MatrixComplexF64,
1049    diag: &mut crate::VectorF64,
1050    subdiag: &mut crate::VectorF64,
1051) -> Result<(), Value> {
1052    let ret = unsafe {
1053        sys::gsl_linalg_hermtd_unpack(
1054            a.unwrap_shared(),
1055            tau.unwrap_shared(),
1056            u.unwrap_unique(),
1057            diag.unwrap_unique(),
1058            subdiag.unwrap_unique(),
1059        )
1060    };
1061    result_handler!(ret, ())
1062}
1063
1064/// This function unpacks the diagonal and subdiagonal of the encoded tridiagonal decomposition (A, tau) obtained from the
1065/// gsl_linalg_hermtd_decomp into the real vectors diag and subdiag.
1066#[doc(alias = "gsl_linalg_hermtd_unpack_T")]
1067pub fn hermtd_unpack_T(
1068    a: &::MatrixComplexF64,
1069    diag: &mut crate::VectorF64,
1070    subdiag: &mut crate::VectorF64,
1071) -> Result<(), Value> {
1072    let ret = unsafe {
1073        sys::gsl_linalg_hermtd_unpack_T(
1074            a.unwrap_shared(),
1075            diag.unwrap_unique(),
1076            subdiag.unwrap_unique(),
1077        )
1078    };
1079    result_handler!(ret, ())
1080}
1081
1082/// This function computes the Hessenberg decomposition of the matrix A by applying the similarity transformation H = U^T A U. On output, H
1083/// 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.
1084/// U is a product of N - 2 Householder matrices. The Householder vectors are stored in the lower portion of A (below the subdiagonal) and
1085/// the Householder coefficients are stored in the vector tau. tau must be of length N.
1086#[doc(alias = "gsl_linalg_hessenberg_decomp")]
1087pub fn hessenberg_decomp(
1088    a: &mut crate::MatrixF64,
1089    tau: &mut crate::VectorF64,
1090) -> Result<(), Value> {
1091    let ret = unsafe { sys::gsl_linalg_hessenberg_decomp(a.unwrap_unique(), tau.unwrap_unique()) };
1092    result_handler!(ret, ())
1093}
1094
1095/// This function constructs the orthogonal matrix U from the information stored in the Hessenberg matrix H along with the vector tau. H and
1096/// tau are outputs from gsl_linalg_hessenberg_decomp.
1097#[doc(alias = "gsl_linalg_hessenberg_unpack")]
1098pub fn hessenberg_unpack(
1099    h: &mut crate::MatrixF64,
1100    tau: &mut crate::VectorF64,
1101    u: &mut crate::MatrixF64,
1102) -> Result<(), Value> {
1103    let ret = unsafe {
1104        sys::gsl_linalg_hessenberg_unpack(h.unwrap_unique(), tau.unwrap_unique(), u.unwrap_unique())
1105    };
1106    result_handler!(ret, ())
1107}
1108
1109/// 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
1110/// be initialized prior to calling this function. Setting V to the identity matrix provides the same result as gsl_linalg_hessenberg_unpack.
1111/// If H is order N, then V must have N columns but may have any number of rows.
1112#[doc(alias = "gsl_linalg_hessenberg_unpack_accum")]
1113pub fn hessenberg_unpack_accum(
1114    h: &mut crate::MatrixF64,
1115    tau: &mut crate::VectorF64,
1116    v: &mut crate::MatrixF64,
1117) -> Result<(), Value> {
1118    let ret = unsafe {
1119        sys::gsl_linalg_hessenberg_unpack_accum(
1120            h.unwrap_unique(),
1121            tau.unwrap_unique(),
1122            v.unwrap_unique(),
1123        )
1124    };
1125    result_handler!(ret, ())
1126}
1127
1128/// This function sets the lower triangular portion of H, below the subdiagonal, to zero. It is useful for clearing out the Householder
1129/// vectors after calling gsl_linalg_hessenberg_decomp.
1130#[doc(alias = "gsl_linalg_hessenberg_set_zero")]
1131pub fn hessenberg_set_zero(h: &mut crate::MatrixF64) -> Result<(), Value> {
1132    let ret = unsafe { sys::gsl_linalg_hessenberg_set_zero(h.unwrap_unique()) };
1133    result_handler!(ret, ())
1134}
1135
1136/// This function computes the Hessenberg-Triangular decomposition of the matrix pair (A, B). On output, H is stored in A, and R is stored
1137/// in B. If U and V are provided (they may be null), the similarity transformations are stored in them. Additional workspace of length N
1138/// is needed in work.
1139#[doc(alias = "gsl_linalg_hesstri_decomp")]
1140pub fn hesstri_decomp(
1141    a: &mut crate::MatrixF64,
1142    b: &mut crate::MatrixF64,
1143    u: &mut crate::MatrixF64,
1144    v: &mut crate::MatrixF64,
1145    work: &mut crate::VectorF64,
1146) -> Result<(), Value> {
1147    let ret = unsafe {
1148        sys::gsl_linalg_hesstri_decomp(
1149            a.unwrap_unique(),
1150            b.unwrap_unique(),
1151            u.unwrap_unique(),
1152            v.unwrap_unique(),
1153            work.unwrap_unique(),
1154        )
1155    };
1156    result_handler!(ret, ())
1157}
1158
1159/// This function factorizes the M-by-N matrix A into bidiagonal form U B V^T. The diagonal and superdiagonal of the matrix B are stored in
1160/// the diagonal and superdiagonal of A. The orthogonal matrices U and V are stored as compressed Householder vectors in the remaining elements
1161/// of A. The Householder coefficients are stored in the vectors tau_U and tau_V. The length of tau_U must equal the number of elements in
1162/// the diagonal of A and the length of tau_V should be one element shorter.
1163#[doc(alias = "gsl_linalg_bidiag_decomp")]
1164pub fn bidiag_decomp(
1165    a: &mut crate::MatrixF64,
1166    tau_u: &mut crate::VectorF64,
1167    tau_v: &mut crate::VectorF64,
1168) -> Result<(), Value> {
1169    let ret = unsafe {
1170        sys::gsl_linalg_bidiag_decomp(
1171            a.unwrap_unique(),
1172            tau_u.unwrap_unique(),
1173            tau_v.unwrap_unique(),
1174        )
1175    };
1176    result_handler!(ret, ())
1177}
1178
1179/// This function unpacks the bidiagonal decomposition of A produced by gsl_linalg_bidiag_decomp, (A, tau_U, tau_V) into the separate orthogonal
1180/// matrices U, V and the diagonal vector diag and superdiagonal superdiag. Note that U is stored as a compact M-by-N orthogonal matrix satisfying
1181/// U^T U = I for efficiency.
1182#[doc(alias = "gsl_linalg_bidiag_unpack")]
1183pub fn bidiag_unpack(
1184    a: &mut crate::MatrixF64,
1185    tau_u: &::VectorF64,
1186    u: &mut crate::MatrixF64,
1187    tau_v: &::VectorF64,
1188    v: &mut crate::MatrixF64,
1189    diag: &mut crate::VectorF64,
1190    superdiag: &mut crate::VectorF64,
1191) -> Result<(), Value> {
1192    let ret = unsafe {
1193        sys::gsl_linalg_bidiag_unpack(
1194            a.unwrap_unique(),
1195            tau_u.unwrap_shared(),
1196            u.unwrap_unique(),
1197            tau_v.unwrap_shared(),
1198            v.unwrap_unique(),
1199            diag.unwrap_unique(),
1200            superdiag.unwrap_unique(),
1201        )
1202    };
1203    result_handler!(ret, ())
1204}
1205
1206/// This function unpacks the bidiagonal decomposition of A produced by gsl_linalg_bidiag_decomp, (A, tau_U, tau_V) into the separate orthogonal
1207/// matrices U, V and the diagonal vector diag and superdiagonal superdiag. The matrix U is stored in-place in A.
1208#[doc(alias = "gsl_linalg_bidiag_unpack2")]
1209pub fn bidiag_unpack2(
1210    a: &mut crate::MatrixF64,
1211    tau_u: &mut crate::VectorF64,
1212    tau_v: &mut crate::VectorF64,
1213    v: &mut crate::MatrixF64,
1214) -> Result<(), Value> {
1215    let ret = unsafe {
1216        sys::gsl_linalg_bidiag_unpack2(
1217            a.unwrap_unique(),
1218            tau_u.unwrap_unique(),
1219            tau_v.unwrap_unique(),
1220            v.unwrap_unique(),
1221        )
1222    };
1223    result_handler!(ret, ())
1224}
1225
1226/// This function unpacks the diagonal and superdiagonal of the bidiagonal decomposition of A from gsl_linalg_bidiag_decomp, into the diagonal
1227/// vector diag and superdiagonal vector superdiag.
1228#[doc(alias = "gsl_linalg_bidiag_unpack_B")]
1229pub fn bidiag_unpack_B(
1230    a: &::MatrixF64,
1231    diag: &mut crate::VectorF64,
1232    superdiag: &mut crate::VectorF64,
1233) -> Result<(), Value> {
1234    let ret = unsafe {
1235        sys::gsl_linalg_bidiag_unpack_B(
1236            a.unwrap_shared(),
1237            diag.unwrap_unique(),
1238            superdiag.unwrap_unique(),
1239        )
1240    };
1241    result_handler!(ret, ())
1242}
1243
1244/// This function prepares a Householder transformation P = I - \tau v v^T which can be used to zero all the elements of the input vector except
1245/// the first. On output the transformation is stored in the vector v and the scalar \tau is returned.
1246#[doc(alias = "gsl_linalg_householder_transform")]
1247pub fn householder_transform(v: &mut crate::VectorF64) -> f64 {
1248    unsafe { sys::gsl_linalg_householder_transform(v.unwrap_unique()) }
1249}
1250
1251/// This function prepares a Householder transformation P = I - \tau v v^T which can be used to zero all the elements of the input vector except
1252/// the first. On output the transformation is stored in the vector v and the scalar \tau is returned.
1253#[doc(alias = "gsl_linalg_complex_householder_transform")]
1254pub fn complex_householder_transform(v: &mut crate::VectorComplexF64) -> crate::ComplexF64 {
1255    unsafe {
1256        std::mem::transmute(sys::gsl_linalg_complex_householder_transform(
1257            v.unwrap_unique(),
1258        ))
1259    }
1260}
1261
1262/// 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
1263/// the result P A is stored in A.
1264#[doc(alias = "gsl_linalg_householder_hm")]
1265pub fn householder_hm(tau: f64, v: &::VectorF64, a: &mut crate::MatrixF64) -> Result<(), Value> {
1266    let ret = unsafe { sys::gsl_linalg_householder_hm(tau, v.unwrap_shared(), a.unwrap_unique()) };
1267    result_handler!(ret, ())
1268}
1269
1270/// 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
1271/// the result P A is stored in A.
1272#[doc(alias = "gsl_linalg_complex_householder_hm")]
1273pub fn complex_householder_hm(
1274    tau: &::ComplexF64,
1275    v: &::VectorComplexF64,
1276    a: &mut crate::MatrixComplexF64,
1277) -> Result<(), Value> {
1278    let ret = unsafe {
1279        sys::gsl_linalg_complex_householder_hm(
1280            std::mem::transmute(*tau),
1281            v.unwrap_shared(),
1282            a.unwrap_unique(),
1283        )
1284    };
1285    result_handler!(ret, ())
1286}
1287
1288/// 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
1289/// the result A P is stored in A.
1290#[doc(alias = "gsl_linalg_householder_mh")]
1291pub fn householder_mh(tau: f64, v: &::VectorF64, a: &mut crate::MatrixF64) -> Result<(), Value> {
1292    let ret = unsafe { sys::gsl_linalg_householder_mh(tau, v.unwrap_shared(), a.unwrap_unique()) };
1293    result_handler!(ret, ())
1294}
1295
1296/// 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
1297/// the result A P is stored in A.
1298#[doc(alias = "gsl_linalg_complex_householder_mh")]
1299pub fn complex_householder_mh(
1300    tau: &::ComplexF64,
1301    v: &::VectorComplexF64,
1302    a: &mut crate::MatrixComplexF64,
1303) -> Result<(), Value> {
1304    let ret = unsafe {
1305        sys::gsl_linalg_complex_householder_mh(
1306            std::mem::transmute(*tau),
1307            v.unwrap_shared(),
1308            a.unwrap_unique(),
1309        )
1310    };
1311    result_handler!(ret, ())
1312}
1313
1314/// 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
1315/// w is stored in w.
1316#[doc(alias = "gsl_linalg_householder_hv")]
1317pub fn householder_hv(tau: f64, v: &::VectorF64, w: &mut crate::VectorF64) -> Result<(), Value> {
1318    let ret = unsafe { sys::gsl_linalg_householder_hv(tau, v.unwrap_shared(), w.unwrap_unique()) };
1319    result_handler!(ret, ())
1320}
1321
1322/// 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
1323/// w is stored in w.
1324#[doc(alias = "gsl_linalg_complex_householder_hv")]
1325pub fn complex_householder_hv(
1326    tau: &::ComplexF64,
1327    v: &::VectorComplexF64,
1328    w: &mut crate::VectorComplexF64,
1329) -> Result<(), Value> {
1330    let ret = unsafe {
1331        sys::gsl_linalg_complex_householder_hv(
1332            std::mem::transmute(*tau),
1333            v.unwrap_shared(),
1334            w.unwrap_unique(),
1335        )
1336    };
1337    result_handler!(ret, ())
1338}
1339
1340/// This function solves the system A x = b directly using Householder transformations. On output the solution is stored in x and b is not
1341/// modified. The matrix A is destroyed by the Householder transformations.
1342#[doc(alias = "gsl_linalg_HH_solve")]
1343pub fn HH_solve(
1344    mut a: crate::MatrixF64,
1345    b: &::VectorF64,
1346    x: &mut crate::VectorF64,
1347) -> Result<(), Value> {
1348    let ret = unsafe {
1349        sys::gsl_linalg_HH_solve(a.unwrap_unique(), b.unwrap_shared(), x.unwrap_unique())
1350    };
1351    result_handler!(ret, ())
1352}
1353
1354/// This function solves the system A x = b in-place using Householder transformations. On input x should contain the right-hand side b,
1355/// which is replaced by the solution on output. The matrix A is destroyed by the Householder transformations.
1356#[doc(alias = "gsl_linalg_HH_svx")]
1357pub fn HH_svx(mut a: crate::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
1358    let ret = unsafe { sys::gsl_linalg_HH_svx(a.unwrap_unique(), x.unwrap_unique()) };
1359    result_handler!(ret, ())
1360}
1361
1362/// This function solves the general N-by-N system A x = b where A is tridiagonal (N >= 2). The super-diagonal and sub-diagonal vectors
1363/// e and f must be one element shorter than the diagonal vector diag. The form of A for the 4-by-4 case is shown below,
1364///
1365/// ```text
1366/// A = ( d_0 e_0  0   0  )
1367///     ( f_0 d_1 e_1  0  )
1368///     (  0  f_1 d_2 e_2 )
1369///     (  0   0  f_2 d_3 )
1370/// ```
1371#[doc(alias = "gsl_linalg_solve_tridiag")]
1372pub fn solve_tridiag(
1373    diag: &::VectorF64,
1374    e: &::VectorF64,
1375    f: &::VectorF64,
1376    b: &::VectorF64,
1377    x: &mut crate::VectorF64,
1378) -> Result<(), Value> {
1379    let ret = unsafe {
1380        sys::gsl_linalg_solve_tridiag(
1381            diag.unwrap_shared(),
1382            e.unwrap_shared(),
1383            f.unwrap_shared(),
1384            b.unwrap_shared(),
1385            x.unwrap_unique(),
1386        )
1387    };
1388    result_handler!(ret, ())
1389}
1390
1391/// This function solves the general N-by-N system A x = b where A is symmetric tridiagonal (N >= 2). The off-diagonal vector e must be one
1392/// element shorter than the diagonal vector diag. The form of A for the 4-by-4 case is shown below,
1393///
1394/// ```text
1395/// A = ( d_0 e_0  0   0  )
1396///     ( e_0 d_1 e_1  0  )
1397///     (  0  e_1 d_2 e_2 )
1398///     (  0   0  e_2 d_3 )
1399/// ```
1400#[doc(alias = "gsl_linalg_solve_symm_tridiag")]
1401pub fn solve_symm_tridiag(
1402    diag: &::VectorF64,
1403    e: &::VectorF64,
1404    b: &::VectorF64,
1405    x: &mut crate::VectorF64,
1406) -> Result<(), Value> {
1407    let ret = unsafe {
1408        sys::gsl_linalg_solve_symm_tridiag(
1409            diag.unwrap_shared(),
1410            e.unwrap_shared(),
1411            b.unwrap_shared(),
1412            x.unwrap_unique(),
1413        )
1414    };
1415    result_handler!(ret, ())
1416}
1417
1418/// This function solves the general N-by-N system A x = b where A is cyclic tridiagonal (N >= 3). The cyclic super-diagonal and sub-diagonal
1419/// vectors e and f must have the same number of elements as the diagonal vector diag. The form of A for the 4-by-4 case is shown below,
1420///
1421/// ```text
1422/// A = ( d_0 e_0  0  f_3 )
1423///     ( f_0 d_1 e_1  0  )
1424///     (  0  f_1 d_2 e_2 )
1425///     ( e_3  0  f_2 d_3 )
1426/// ```
1427#[doc(alias = "gsl_linalg_solve_cyc_tridiag")]
1428pub fn solve_cyc_tridiag(
1429    diag: &::VectorF64,
1430    e: &::VectorF64,
1431    f: &::VectorF64,
1432    b: &::VectorF64,
1433    x: &mut crate::VectorF64,
1434) -> Result<(), Value> {
1435    let ret = unsafe {
1436        sys::gsl_linalg_solve_cyc_tridiag(
1437            diag.unwrap_shared(),
1438            e.unwrap_shared(),
1439            f.unwrap_shared(),
1440            b.unwrap_shared(),
1441            x.unwrap_unique(),
1442        )
1443    };
1444    result_handler!(ret, ())
1445}
1446
1447/// This function solves the general N-by-N system A x = b where A is symmetric cyclic tridiagonal (N >= 3). The cyclic off-diagonal vector
1448/// e must have the same number of elements as the diagonal vector diag. The form of A for the 4-by-4 case is shown below,
1449///
1450/// ```text
1451/// A = ( d_0 e_0  0  e_3 )
1452///     ( e_0 d_1 e_1  0  )
1453///     (  0  e_1 d_2 e_2 )
1454///     ( e_3  0  e_2 d_3 )
1455/// ```
1456#[doc(alias = "gsl_linalg_solve_symm_cyc_tridiag")]
1457pub fn solve_symm_cyc_tridiag(
1458    diag: &::VectorF64,
1459    e: &::VectorF64,
1460    b: &::VectorF64,
1461    x: &mut crate::VectorF64,
1462) -> Result<(), Value> {
1463    let ret = unsafe {
1464        sys::gsl_linalg_solve_symm_cyc_tridiag(
1465            diag.unwrap_shared(),
1466            e.unwrap_shared(),
1467            b.unwrap_shared(),
1468            x.unwrap_unique(),
1469        )
1470    };
1471    result_handler!(ret, ())
1472}
1473
1474/// This function replaces the matrix A with its balanced counterpart and stores the diagonal elements of the similarity transformation into
1475/// the vector D.
1476#[doc(alias = "gsl_linalg_balance_matrix")]
1477pub fn balance_matrix(a: &mut crate::MatrixF64, d: &mut crate::VectorF64) -> Result<(), Value> {
1478    let ret = unsafe { sys::gsl_linalg_balance_matrix(a.unwrap_unique(), d.unwrap_unique()) };
1479    result_handler!(ret, ())
1480}
1481
1482#[cfg(feature = "v2_2")]
1483#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1484#[doc(alias = "gsl_linalg_pcholesky_decomp")]
1485pub fn pcholesky_decomp(a: &mut crate::MatrixF64, p: &mut ::Permutation) -> Result<(), Value> {
1486    let ret = unsafe { sys::gsl_linalg_pcholesky_decomp(a.unwrap_unique(), p.unwrap_unique()) };
1487    result_handler!(ret, ())
1488}
1489
1490#[cfg(feature = "v2_2")]
1491#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1492#[doc(alias = "gsl_linalg_pcholesky_solve")]
1493pub fn pcholesky_solve(
1494    LDLT: &::MatrixF64,
1495    p: &::Permutation,
1496    b: &::VectorF64,
1497    x: &mut crate::VectorF64,
1498) -> Result<(), Value> {
1499    let ret = unsafe {
1500        sys::gsl_linalg_pcholesky_solve(
1501            LDLT.unwrap_shared(),
1502            p.unwrap_shared(),
1503            b.unwrap_shared(),
1504            x.unwrap_unique(),
1505        )
1506    };
1507    result_handler!(ret, ())
1508}
1509
1510#[cfg(feature = "v2_2")]
1511#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1512#[doc(alias = "gsl_linalg_pcholesky_svx")]
1513pub fn pcholesky_svx(
1514    LDLT: &::MatrixF64,
1515    p: &::Permutation,
1516    x: &mut crate::VectorF64,
1517) -> Result<(), Value> {
1518    let ret = unsafe {
1519        sys::gsl_linalg_pcholesky_svx(LDLT.unwrap_shared(), p.unwrap_shared(), x.unwrap_unique())
1520    };
1521    result_handler!(ret, ())
1522}
1523
1524#[cfg(feature = "v2_2")]
1525#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1526#[doc(alias = "gsl_linalg_pcholesky_decomp2")]
1527pub fn pcholesky_decomp2(
1528    A: &mut crate::MatrixF64,
1529    p: &mut ::Permutation,
1530    S: &mut crate::VectorF64,
1531) -> Result<(), Value> {
1532    let ret = unsafe {
1533        sys::gsl_linalg_pcholesky_decomp2(A.unwrap_unique(), p.unwrap_unique(), S.unwrap_unique())
1534    };
1535    result_handler!(ret, ())
1536}
1537
1538#[cfg(feature = "v2_2")]
1539#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1540#[doc(alias = "gsl_linalg_pcholesky_solve2")]
1541pub fn pcholesky_solve2(
1542    LDLT: &::MatrixF64,
1543    p: &::Permutation,
1544    S: &::VectorF64,
1545    b: &::VectorF64,
1546    x: &mut crate::VectorF64,
1547) -> Result<(), Value> {
1548    let ret = unsafe {
1549        sys::gsl_linalg_pcholesky_solve2(
1550            LDLT.unwrap_shared(),
1551            p.unwrap_shared(),
1552            S.unwrap_shared(),
1553            b.unwrap_shared(),
1554            x.unwrap_unique(),
1555        )
1556    };
1557    result_handler!(ret, ())
1558}
1559
1560#[cfg(feature = "v2_2")]
1561#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1562#[doc(alias = "gsl_linalg_pcholesky_svx2")]
1563pub fn pcholesky_svx2(
1564    LDLT: &::MatrixF64,
1565    p: &::Permutation,
1566    S: &::VectorF64,
1567    x: &mut crate::VectorF64,
1568) -> Result<(), Value> {
1569    let ret = unsafe {
1570        sys::gsl_linalg_pcholesky_svx2(
1571            LDLT.unwrap_shared(),
1572            p.unwrap_shared(),
1573            S.unwrap_shared(),
1574            x.unwrap_unique(),
1575        )
1576    };
1577    result_handler!(ret, ())
1578}
1579
1580#[cfg(feature = "v2_2")]
1581#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1582#[doc(alias = "gsl_linalg_pcholesky_invert")]
1583pub fn pcholesky_invert(
1584    LDLT: &::MatrixF64,
1585    p: &::Permutation,
1586    Ainv: &mut crate::MatrixF64,
1587) -> Result<(), Value> {
1588    let ret = unsafe {
1589        sys::gsl_linalg_pcholesky_invert(
1590            LDLT.unwrap_shared(),
1591            p.unwrap_shared(),
1592            Ainv.unwrap_unique(),
1593        )
1594    };
1595    result_handler!(ret, ())
1596}
1597
1598/// Returns `(Value, rcond)`.
1599#[cfg(feature = "v2_2")]
1600#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1601#[doc(alias = "gsl_linalg_pcholesky_rcond")]
1602pub fn pcholesky_rcond(
1603    LDLT: &::MatrixF64,
1604    p: &::Permutation,
1605    work: &mut crate::VectorF64,
1606) -> Result<f64, Value> {
1607    let mut rcond = 0.;
1608    let ret = unsafe {
1609        sys::gsl_linalg_pcholesky_rcond(
1610            LDLT.unwrap_shared(),
1611            p.unwrap_shared(),
1612            &mut rcond,
1613            work.unwrap_unique(),
1614        )
1615    };
1616    result_handler!(ret, rcond)
1617}
1618
1619#[cfg(feature = "v2_2")]
1620#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1621#[doc(alias = "gsl_linalg_mcholesky_decomp")]
1622pub fn mcholesky_decomp(
1623    A: &mut crate::MatrixF64,
1624    p: &mut ::Permutation,
1625    E: &mut crate::VectorF64,
1626) -> Result<(), Value> {
1627    let ret = unsafe {
1628        sys::gsl_linalg_mcholesky_decomp(A.unwrap_unique(), p.unwrap_unique(), E.unwrap_unique())
1629    };
1630    result_handler!(ret, ())
1631}
1632
1633#[cfg(feature = "v2_2")]
1634#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1635#[doc(alias = "gsl_linalg_mcholesky_solve")]
1636pub fn mcholesky_solve(
1637    LDLT: &::MatrixF64,
1638    p: &::Permutation,
1639    b: &::VectorF64,
1640    x: &mut crate::VectorF64,
1641) -> Result<(), Value> {
1642    let ret = unsafe {
1643        sys::gsl_linalg_mcholesky_solve(
1644            LDLT.unwrap_shared(),
1645            p.unwrap_shared(),
1646            b.unwrap_shared(),
1647            x.unwrap_unique(),
1648        )
1649    };
1650    result_handler!(ret, ())
1651}
1652
1653#[cfg(feature = "v2_2")]
1654#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1655#[doc(alias = "gsl_linalg_mcholesky_svx")]
1656pub fn mcholesky_svx(
1657    LDLT: &::MatrixF64,
1658    p: &::Permutation,
1659    x: &mut crate::VectorF64,
1660) -> Result<(), Value> {
1661    let ret = unsafe {
1662        sys::gsl_linalg_mcholesky_svx(LDLT.unwrap_shared(), p.unwrap_shared(), x.unwrap_unique())
1663    };
1664    result_handler!(ret, ())
1665}
1666
1667/// Returns `(Value, rcond)`.
1668#[cfg(feature = "v2_2")]
1669#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1670#[doc(alias = "gsl_linalg_mcholesky_rcond")]
1671pub fn mcholesky_rcond(
1672    LDLT: &::MatrixF64,
1673    p: &::Permutation,
1674    work: &mut crate::VectorF64,
1675) -> Result<f64, Value> {
1676    let mut rcond = 0.;
1677    let ret = unsafe {
1678        sys::gsl_linalg_mcholesky_rcond(
1679            LDLT.unwrap_shared(),
1680            p.unwrap_shared(),
1681            &mut rcond,
1682            work.unwrap_unique(),
1683        )
1684    };
1685    result_handler!(ret, rcond)
1686}
1687
1688#[cfg(feature = "v2_2")]
1689#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1690#[doc(alias = "gsl_linalg_mcholesky_invert")]
1691pub fn mcholesky_invert(
1692    LDLT: &::MatrixF64,
1693    p: &::Permutation,
1694    Ainv: &mut crate::MatrixF64,
1695) -> Result<(), Value> {
1696    let ret = unsafe {
1697        sys::gsl_linalg_mcholesky_invert(
1698            LDLT.unwrap_shared(),
1699            p.unwrap_shared(),
1700            Ainv.unwrap_unique(),
1701        )
1702    };
1703    result_handler!(ret, ())
1704}
1705
1706#[cfg(feature = "v2_6")]
1707#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1708#[doc(alias = "gsl_linalg_cholesky_band_decomp")]
1709pub fn cholesky_band_decomp(A: &mut crate::MatrixF64) -> Result<(), Value> {
1710    let ret = unsafe { sys::gsl_linalg_cholesky_band_decomp(A.unwrap_unique()) };
1711    result_handler!(ret, ())
1712}
1713
1714#[cfg(feature = "v2_6")]
1715#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1716#[doc(alias = "gsl_linalg_cholesky_band_solve")]
1717pub fn cholesky_band_solve(
1718    LLT: &::MatrixF64,
1719    b: &::VectorF64,
1720    x: &mut crate::VectorF64,
1721) -> Result<(), Value> {
1722    let ret = unsafe {
1723        sys::gsl_linalg_cholesky_band_solve(
1724            LLT.unwrap_shared(),
1725            b.unwrap_shared(),
1726            x.unwrap_unique(),
1727        )
1728    };
1729    result_handler!(ret, ())
1730}
1731
1732#[cfg(feature = "v2_6")]
1733#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1734#[doc(alias = "gsl_linalg_cholesky_band_svx")]
1735pub fn cholesky_band_svx(LLT: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
1736    let ret = unsafe { sys::gsl_linalg_cholesky_band_svx(LLT.unwrap_shared(), x.unwrap_unique()) };
1737    result_handler!(ret, ())
1738}
1739
1740#[cfg(feature = "v2_7")]
1741#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_7")))]
1742#[doc(alias = "gsl_linalg_cholesky_band_solvem")]
1743pub fn cholesky_band_solvem(
1744    LLT: &::MatrixF64,
1745    B: &::MatrixF64,
1746    X: &mut crate::MatrixF64,
1747) -> Result<(), Value> {
1748    let ret = unsafe {
1749        sys::gsl_linalg_cholesky_band_solvem(
1750            LLT.unwrap_shared(),
1751            B.unwrap_shared(),
1752            X.unwrap_unique(),
1753        )
1754    };
1755    result_handler!(ret, ())
1756}
1757
1758#[cfg(feature = "v2_7")]
1759#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_7")))]
1760#[doc(alias = "gsl_linalg_cholesky_band_svxm")]
1761pub fn cholesky_band_svxm(LLT: &::MatrixF64, X: &mut crate::MatrixF64) -> Result<(), Value> {
1762    let ret = unsafe { sys::gsl_linalg_cholesky_band_svxm(LLT.unwrap_shared(), X.unwrap_unique()) };
1763    result_handler!(ret, ())
1764}
1765
1766#[cfg(feature = "v2_6")]
1767#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1768#[doc(alias = "gsl_linalg_cholesky_band_invert")]
1769pub fn cholesky_band_invert(LLT: &::MatrixF64, Ainv: &mut crate::MatrixF64) -> Result<(), Value> {
1770    let ret =
1771        unsafe { sys::gsl_linalg_cholesky_band_invert(LLT.unwrap_shared(), Ainv.unwrap_unique()) };
1772    result_handler!(ret, ())
1773}
1774
1775#[cfg(feature = "v2_6")]
1776#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1777#[doc(alias = "gsl_linalg_cholesky_band_unpack")]
1778pub fn cholesky_band_unpack(LLT: &::MatrixF64, L: &mut crate::MatrixF64) -> Result<(), Value> {
1779    let ret =
1780        unsafe { sys::gsl_linalg_cholesky_band_unpack(LLT.unwrap_shared(), L.unwrap_unique()) };
1781    result_handler!(ret, ())
1782}
1783
1784/// Returns `(Value, rcond)`.
1785#[cfg(feature = "v2_6")]
1786#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1787#[doc(alias = "gsl_linalg_cholesky_band_rcond")]
1788pub fn cholesky_band_rcond(LLT: &::MatrixF64, work: &mut crate::VectorF64) -> Result<f64, Value> {
1789    let mut rcond = 0.;
1790    let ret = unsafe {
1791        sys::gsl_linalg_cholesky_band_rcond(LLT.unwrap_shared(), &mut rcond, work.unwrap_unique())
1792    };
1793    result_handler!(ret, rcond)
1794}
1795
1796#[cfg(feature = "v2_6")]
1797#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1798#[doc(alias = "gsl_linalg_ldlt_decomp")]
1799pub fn ldlt_decomp(A: &mut crate::MatrixF64) -> Result<(), Value> {
1800    let ret = unsafe { sys::gsl_linalg_ldlt_decomp(A.unwrap_unique()) };
1801    result_handler!(ret, ())
1802}
1803
1804#[cfg(feature = "v2_6")]
1805#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1806#[doc(alias = "gsl_linalg_ldlt_solve")]
1807pub fn ldlt_solve(
1808    LDLT: &::MatrixF64,
1809    b: &::VectorF64,
1810    x: &mut crate::VectorF64,
1811) -> Result<(), Value> {
1812    let ret = unsafe {
1813        sys::gsl_linalg_ldlt_solve(LDLT.unwrap_shared(), b.unwrap_shared(), x.unwrap_unique())
1814    };
1815    result_handler!(ret, ())
1816}
1817
1818#[cfg(feature = "v2_6")]
1819#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1820#[doc(alias = "gsl_linalg_ldlt_svx")]
1821pub fn ldlt_svx(LDLT: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
1822    let ret = unsafe { sys::gsl_linalg_ldlt_svx(LDLT.unwrap_shared(), x.unwrap_unique()) };
1823    result_handler!(ret, ())
1824}
1825
1826/// Returns `(Value, rcond)`.
1827#[cfg(feature = "v2_6")]
1828#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1829#[doc(alias = "gsl_linalg_ldlt_rcond")]
1830pub fn ldlt_rcond(LDLT: &::MatrixF64, work: &mut crate::VectorF64) -> Result<f64, Value> {
1831    let mut rcond = 0.;
1832    let ret = unsafe {
1833        sys::gsl_linalg_ldlt_rcond(LDLT.unwrap_shared(), &mut rcond, work.unwrap_unique())
1834    };
1835    result_handler!(ret, rcond)
1836}
1837
1838#[cfg(feature = "v2_6")]
1839#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1840#[doc(alias = "gsl_linalg_ldlt_band_decomp")]
1841pub fn ldlt_band_decomp(A: &mut crate::MatrixF64) -> Result<(), Value> {
1842    let ret = unsafe { sys::gsl_linalg_ldlt_band_decomp(A.unwrap_unique()) };
1843    result_handler!(ret, ())
1844}
1845
1846#[cfg(feature = "v2_6")]
1847#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1848#[doc(alias = "gsl_linalg_ldlt_band_solve")]
1849pub fn ldlt_band_solve(
1850    LDLT: &::MatrixF64,
1851    b: &::VectorF64,
1852    x: &mut crate::VectorF64,
1853) -> Result<(), Value> {
1854    let ret = unsafe {
1855        sys::gsl_linalg_ldlt_band_solve(LDLT.unwrap_shared(), b.unwrap_shared(), x.unwrap_unique())
1856    };
1857    result_handler!(ret, ())
1858}
1859
1860#[cfg(feature = "v2_6")]
1861#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1862#[doc(alias = "gsl_linalg_ldlt_band_svx")]
1863pub fn ldlt_band_svx(LDLT: &::MatrixF64, x: &mut crate::VectorF64) -> Result<(), Value> {
1864    let ret = unsafe { sys::gsl_linalg_ldlt_band_svx(LDLT.unwrap_shared(), x.unwrap_unique()) };
1865    result_handler!(ret, ())
1866}
1867
1868#[cfg(feature = "v2_6")]
1869#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1870#[doc(alias = "gsl_linalg_ldlt_band_unpack")]
1871pub fn ldlt_band_unpack(
1872    LDLT: &::MatrixF64,
1873    L: &mut crate::MatrixF64,
1874    D: &mut crate::VectorF64,
1875) -> Result<(), Value> {
1876    let ret = unsafe {
1877        sys::gsl_linalg_ldlt_band_unpack(LDLT.unwrap_shared(), L.unwrap_unique(), D.unwrap_unique())
1878    };
1879    result_handler!(ret, ())
1880}
1881
1882/// Returns `(Value, rcond)`.
1883#[cfg(feature = "v2_6")]
1884#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_6")))]
1885#[doc(alias = "gsl_linalg_ldlt_band_rcond")]
1886pub fn ldlt_band_rcond(LDLT: &::MatrixF64, work: &mut crate::VectorF64) -> Result<f64, Value> {
1887    let mut rcond = 0.;
1888    let ret = unsafe {
1889        sys::gsl_linalg_ldlt_band_rcond(LDLT.unwrap_shared(), &mut rcond, work.unwrap_unique())
1890    };
1891    result_handler!(ret, rcond)
1892}
1893
1894#[cfg(feature = "v2_2")]
1895#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1896#[doc(alias = "gsl_linalg_tri_upper_invert")]
1897pub fn tri_upper_invert(T: &mut crate::MatrixF64) -> Result<(), Value> {
1898    let ret = unsafe { sys::gsl_linalg_tri_upper_invert(T.unwrap_unique()) };
1899    result_handler!(ret, ())
1900}
1901
1902#[cfg(feature = "v2_2")]
1903#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1904#[doc(alias = "gsl_linalg_tri_lower_invert")]
1905pub fn tri_lower_invert(T: &mut crate::MatrixF64) -> Result<(), Value> {
1906    let ret = unsafe { sys::gsl_linalg_tri_lower_invert(T.unwrap_unique()) };
1907    result_handler!(ret, ())
1908}
1909
1910#[cfg(feature = "v2_2")]
1911#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1912#[doc(alias = "gsl_linalg_tri_upper_unit_invert")]
1913pub fn tri_upper_unit_invert(T: &mut crate::MatrixF64) -> Result<(), Value> {
1914    let ret = unsafe { sys::gsl_linalg_tri_upper_unit_invert(T.unwrap_unique()) };
1915    result_handler!(ret, ())
1916}
1917
1918#[cfg(feature = "v2_2")]
1919#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1920#[doc(alias = "gsl_linalg_tri_lower_unit_invert")]
1921pub fn tri_lower_unit_invert(T: &mut crate::MatrixF64) -> Result<(), Value> {
1922    let ret = unsafe { sys::gsl_linalg_tri_lower_unit_invert(T.unwrap_unique()) };
1923    result_handler!(ret, ())
1924}
1925
1926#[cfg(feature = "v2_2")]
1927#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1928#[doc(alias = "gsl_linalg_complex_tri_invert")]
1929pub fn tri_invert(
1930    Uplo: enums::CblasUplo,
1931    Diag: enums::CblasDiag,
1932    T: &mut crate::MatrixComplexF64,
1933) -> Result<(), Value> {
1934    let ret =
1935        unsafe { sys::gsl_linalg_complex_tri_invert(Uplo.into(), Diag.into(), T.unwrap_unique()) };
1936    result_handler!(ret, ())
1937}
1938
1939#[doc(alias = "gsl_linalg_complex_tri_invert")]
1940pub fn complex_tri_invert(
1941    Uplo: enums::CblasUplo,
1942    Diag: enums::CblasDiag,
1943    T: &mut crate::MatrixComplexF64,
1944) -> Result<(), Value> {
1945    let ret =
1946        unsafe { sys::gsl_linalg_complex_tri_invert(Uplo.into(), Diag.into(), T.unwrap_unique()) };
1947    result_handler!(ret, ())
1948}
1949
1950#[cfg(feature = "v2_2")]
1951#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1952#[doc(alias = "gsl_linalg_tri_LTL")]
1953pub fn tri_LTL(L: &mut crate::MatrixF64) -> Result<(), Value> {
1954    let ret = unsafe { sys::gsl_linalg_tri_LTL(L.unwrap_unique()) };
1955    result_handler!(ret, ())
1956}
1957
1958#[cfg(feature = "v2_2")]
1959#[cfg_attr(feature = "dox", doc(cfg(feature = "v2_2")))]
1960#[doc(alias = "gsl_linalg_tri_UL")]
1961pub fn tri_UL(LU: &mut crate::MatrixF64) -> Result<(), Value> {
1962    let ret = unsafe { sys::gsl_linalg_tri_UL(LU.unwrap_unique()) };
1963    result_handler!(ret, ())
1964}
1965
1966#[doc(alias = "gsl_linalg_complex_tri_LHL")]
1967pub fn complex_tri_LHL(L: &mut crate::MatrixComplexF64) -> Result<(), Value> {
1968    let ret = unsafe { sys::gsl_linalg_complex_tri_LHL(L.unwrap_unique()) };
1969    result_handler!(ret, ())
1970}
1971
1972#[doc(alias = "gsl_linalg_complex_tri_UL")]
1973pub fn complex_tri_UL(LU: &mut crate::MatrixComplexF64) -> Result<(), Value> {
1974    let ret = unsafe { sys::gsl_linalg_complex_tri_UL(LU.unwrap_unique()) };
1975    result_handler!(ret, ())
1976}
1977
1978/// Returns `(c, s)`.
1979#[doc(alias = "gsl_linalg_givens")]
1980pub fn givens(a: f64, b: f64) -> (f64, f64) {
1981    let mut c = 0.;
1982    let mut s = 0.;
1983    unsafe { sys::gsl_linalg_givens(a, b, &mut c, &mut s) };
1984    (c, s)
1985}
1986
1987#[doc(alias = "gsl_linalg_givens_gv")]
1988pub fn givens_gv(v: &mut crate::VectorF64, i: usize, j: usize, c: f64, s: f64) {
1989    unsafe { sys::gsl_linalg_givens_gv(v.unwrap_unique(), i, j, c, s) }
1990}