20#ifndef CCGSL_LINALG_HPP
21#define CCGSL_LINALG_HPP
24#include<gsl/gsl_linalg.h>
44 return gsl_linalg_householder_transform( v.
get() ); }
52 return gsl_linalg_householder_transform2( &alpha, v.
get() ); }
59 return gsl_linalg_complex_householder_transform( v.
get() ); }
68 return gsl_linalg_householder_hm( tau, v.
get(), A.
get() ); }
77 return gsl_linalg_householder_mh( tau, v.
get(), A.
get() ); }
86 return gsl_linalg_householder_hv( tau, v.
get(), w.
get() ); }
94 return gsl_linalg_householder_hm1( tau, A.
get() ); }
103 return gsl_linalg_complex_householder_hm( tau.
get(), v.
get(), A.
get() ); }
112 return gsl_linalg_complex_householder_mh( tau.
get(), v.
get(), A.
get() ); }
121 return gsl_linalg_complex_householder_hv( tau.
get(), v.
get(), w.
get() ); }
129 return gsl_linalg_hessenberg_decomp( A.
get(), tau.
get() ); }
138 return gsl_linalg_hessenberg_unpack( H.
get(), tau.
get(), U.
get() ); }
147 return gsl_linalg_hessenberg_unpack_accum( H.
get(), tau.
get(), U.
get() ); }
163 return gsl_linalg_hessenberg_submatrix( M.
get(), A.
get(), top, tau.
get() ); }
183 return gsl_linalg_hesstri_decomp( A.
get(), B.
get(), U.
get(), V.
get(), work.
get() ); }
193 return gsl_linalg_SV_decomp( A.
get(), V.
get(), S.
get(), work.
get() ); }
204 return gsl_linalg_SV_decomp_mod( A.
get(), X.
get(), V.
get(), S.
get(), work.
get() ); }
213 return gsl_linalg_SV_decomp_jacobi( A.
get(), Q.
get(), S.
get() ); }
224 return gsl_linalg_SV_solve( U.
get(), Q.
get(), S.
get(),
b.get(), x.
get() ); }
234 return gsl_linalg_LU_decomp( A.
get(), p.
get(), signum ); }
244 return gsl_linalg_LU_decomp( A.
get(), p.
get(), &signum ); }
254 return gsl_linalg_LU_solve( LU.
get(), p.
get(),
b.get(), x.
get() ); }
263 return gsl_linalg_LU_svx( LU.
get(), p.
get(), x.
get() ); }
292 inline double LU_det(
matrix& LU,
int signum ){
return gsl_linalg_LU_det( LU.
get(), signum ); }
316 return gsl_linalg_complex_LU_decomp( A.
get(), p.
get(), signum ); }
326 return gsl_linalg_complex_LU_decomp( A.
get(), p.
get(), &signum ); }
337 return gsl_linalg_complex_LU_solve( LU.
get(), p.
get(),
b.get(), x.
get() ); }
346 return gsl_linalg_complex_LU_svx( LU.
get(), p.
get(), x.
get() ); }
378 return gsl_linalg_complex_LU_det( LU.
get(), signum ); }
385 return gsl_linalg_complex_LU_lndet( LU.
get() ); }
393 return gsl_linalg_complex_LU_sgndet( LU.
get(), signum ); }
410 return gsl_linalg_QR_solve( QR.
get(), tau.
get(),
b.get(), x.
get() ); }
419 return gsl_linalg_QR_svx( QR.
get(), tau.
get(), x.
get() ); }
441 return gsl_linalg_QR_QRsolve( Q.
get(), R.
get(),
b.get(), x.
get() ); }
450 return gsl_linalg_QR_Rsolve( QR.
get(),
b.get(), x.
get() ); }
467 return gsl_linalg_QR_update( Q.
get(), R.
get(), w.
get(), v.
get() ); }
476 return gsl_linalg_QR_QTvec( QR.
get(), tau.
get(), v.
get() ); }
485 return gsl_linalg_QR_Qvec( QR.
get(), tau.
get(), v.
get() ); }
494 return gsl_linalg_QR_QTmat( QR.
get(), tau.
get(), A.
get() ); }
504 return gsl_linalg_QR_unpack( QR.
get(), tau.
get(), Q.
get(), R.
get() ); }
513 return gsl_linalg_R_solve( R.
get(),
b.get(), x.
get() ); }
521 return gsl_linalg_R_svx( R.
get(), x.
get() ); }
533 return gsl_linalg_QRPT_decomp( A.
get(), tau.
get(), p.
get(), signum, norm.
get() ); }
545 return gsl_linalg_QRPT_decomp( A.
get(), tau.
get(), p.
get(), &signum, norm.
get() ); }
560 return gsl_linalg_QRPT_decomp2( A.
get(), q.
get(), r.
get(), tau.
get(), p.
get(), signum, norm.
get() ); }
575 return gsl_linalg_QRPT_decomp2( A.
get(), q.
get(), r.
get(), tau.
get(), p.
get(), &signum, norm.
get() ); }
587 return gsl_linalg_QRPT_solve( QR.
get(), tau.
get(), p.
get(),
b.get(), x.
get() ); }
597 return gsl_linalg_QRPT_svx( QR.
get(), tau.
get(), p.
get(), x.
get() ); }
609 return gsl_linalg_QRPT_QRsolve( Q.
get(), R.
get(), p.
get(),
b.get(), x.
get() ); }
619 return gsl_linalg_QRPT_Rsolve( QR.
get(), p.
get(),
b.get(), x.
get() ); }
628 return gsl_linalg_QRPT_Rsvx( QR.
get(), p.
get(), x.
get() ); }
639 return gsl_linalg_QRPT_update( Q.
get(), R.
get(), p.
get(), u.
get(), v.
get() ); }
647 return gsl_linalg_LQ_decomp( A.
get(), tau.
get() ); }
657 return gsl_linalg_LQ_solve_T( LQ.
get(), tau.
get(),
b.get(), x.
get() ); }
666 return gsl_linalg_LQ_svx_T( LQ.
get(), tau.
get(), x.
get() ); }
687 return gsl_linalg_LQ_Lsolve_T( LQ.
get(),
b.get(), x.
get() ); }
695 return gsl_linalg_LQ_Lsvx_T( LQ.
get(), x.
get() ); }
704 return gsl_linalg_L_solve_T( L.
get(),
b.get(), x.
get() ); }
713 return gsl_linalg_LQ_vecQ( LQ.
get(), tau.
get(), v.
get() ); }
722 return gsl_linalg_LQ_vecQT( LQ.
get(), tau.
get(), v.
get() ); }
732 return gsl_linalg_LQ_unpack( LQ.
get(), tau.
get(), Q.
get(), L.
get() ); }
742 return gsl_linalg_LQ_update( Q.
get(), R.
get(), v.
get(), w.
get() ); }
752 return gsl_linalg_LQ_LQsolve( Q.
get(), L.
get(),
b.get(), x.
get() ); }
764 return gsl_linalg_PTLQ_decomp( A.
get(), tau.
get(), p.
get(), signum, norm.
get() ); }
776 return gsl_linalg_PTLQ_decomp( A.
get(), tau.
get(), p.
get(), &signum, norm.
get() ); }
790 return gsl_linalg_PTLQ_decomp2( A.
get(), q.
get(), r.
get(), tau.
get(), p.
get(), signum, norm.
get() ); }
804 return gsl_linalg_PTLQ_decomp2( A.
get(), q.
get(), r.
get(), tau.
get(), p.
get(), &signum, norm.
get() ); }
816 return gsl_linalg_PTLQ_solve_T( QR.
get(), tau.
get(), p.
get(),
b.get(), x.
get() ); }
826 return gsl_linalg_PTLQ_svx_T( LQ.
get(), tau.
get(), p.
get(), x.
get() ); }
838 return gsl_linalg_PTLQ_LQsolve_T( Q.
get(), L.
get(), p.
get(),
b.get(), x.
get() ); }
848 return gsl_linalg_PTLQ_Lsolve_T( LQ.
get(), p.
get(),
b.get(), x.
get() ); }
857 return gsl_linalg_PTLQ_Lsvx_T( LQ.
get(), p.
get(), x.
get() ); }
868 return gsl_linalg_PTLQ_update( Q.
get(), L.
get(), p.
get(), v.
get(), w.
get() ); }
883 return gsl_linalg_cholesky_solve( cholesky.
get(),
b.get(), x.
get() ); }
891 return gsl_linalg_cholesky_svx( cholesky.
get(), x.
get() ); }
905 return gsl_linalg_cholesky_decomp_unit( A.
get(),
D.get() ); }
912 return gsl_linalg_complex_cholesky_decomp( A.
get() ); }
922 return gsl_linalg_complex_cholesky_solve( cholesky.
get(),
b.get(), x.
get() ); }
930 return gsl_linalg_complex_cholesky_svx( cholesky.
get(), x.
get() ); }
938 return gsl_linalg_symmtd_decomp( A.
get(), tau.
get() ); }
950 return gsl_linalg_symmtd_unpack( A.
get(), tau.
get(), Q.
get(), diag.
get(), subdiag.
get() ); }
959 return gsl_linalg_symmtd_unpack_T( A.
get(), diag.
get(), subdiag.
get() ); }
967 return gsl_linalg_hermtd_decomp( A.
get(), tau.
get() ); }
979 return gsl_linalg_hermtd_unpack( A.
get(), tau.
get(), U.
get(), diag.
get(), sudiag.
get() ); }
988 return gsl_linalg_hermtd_unpack_T( A.
get(), diag.
get(), subdiag.
get() ); }
997 return gsl_linalg_HH_solve( A.
get(),
b.get(), x.
get() ); }
1015 return gsl_linalg_solve_symm_tridiag( diag.
get(), offdiag.
get(),
b.get(), x.
get() ); }
1027 return gsl_linalg_solve_tridiag( diag.
get(), abovediag.
get(), belowdiag.
get(),
b.get(), x.
get() ); }
1038 return gsl_linalg_solve_symm_cyc_tridiag( diag.
get(), offdiag.
get(),
b.get(), x.
get() ); }
1050 return gsl_linalg_solve_cyc_tridiag( diag.
get(), abovediag.
get(), belowdiag.
get(),
b.get(), x.
get() ); }
1059 return gsl_linalg_bidiag_decomp( A.
get(), tau_U.
get(), tau_V.
get() ); }
1073 return gsl_linalg_bidiag_unpack( A.
get(), tau_U.
get(), U.
get(), tau_V.
get(),
1074 V.
get(), diag.
get(), superdiag.
get() ); }
1084 return gsl_linalg_bidiag_unpack2( A.
get(), tau_U.
get(), tau_V.
get(), V.
get() ); }
1093 return gsl_linalg_bidiag_unpack_B( A.
get(), diag.
get(), superdiag.
get() ); }
1101 return gsl_linalg_balance_matrix( A.
get(),
D.get() ); }
1109 return gsl_linalg_balance_accum( A.
get(),
D.get() ); }
1117 return gsl_linalg_balance_columns( A.
get(),
D.get() ); }
1125 return gsl_linalg_SV_leverage( U.
get(), h.
get() ); }
1136 return gsl_linalg_householder_left( tau, v.
get(), A.
get(), work.
get() ); }
1146 return gsl_linalg_householder_right( tau, v.
get(), A.
get(), work.
get() ); }
1157 return gsl_linalg_complex_householder_left( tau.
get(), v.
get(), A.
get(), work.
get() ); }
1165 return gsl_linalg_LU_invx( LU.
get(), p.
get() ); }
1173 return gsl_linalg_complex_LU_invx( LU.
get(), p.
get() ); }
1183 return gsl_linalg_QR_QTmat_r( QR.
get(), T.
get(), B.
get(), work.
get() ); }
1192 return gsl_linalg_QR_matQ( QR.
get(), tau.
get(), A.
get() ); }
1202 return gsl_linalg_QR_unpack_r( QR.
get(), T.
get(), Q.
get(), R.
get() ); }
1211 return gsl_linalg_QR_rcond( QR.
get(), &
rcond, work.
get() ); }
1224 return gsl_linalg_QRPT_lssolve( QR.
get(), tau.
get(), p.
get(),
b.get(), x.
get(),
1239 return gsl_linalg_QRPT_lssolve2( QR.
get(), tau.
get(), p.
get(),
b.get(), rank, x.
get(),
1248 return gsl_linalg_QRPT_rank( QR.
get(), tol ); }
1257 return gsl_linalg_QRPT_rcond( QR.
get(), &
rcond, work.
get() ); }
1269 size_t& rank,
vector& work ){
1270 return gsl_linalg_COD_decomp( A.
get(), tau_Q.
get(), tau_Z.
get(), p.
get(), &rank,
1284 double tol,
size_t& rank,
vector& work ){
1285 return gsl_linalg_COD_decomp_e( A.
get(), tau_Q.
get(), tau_Z.
get(), p.
get(), tol,
1286 &rank, work.
get() ); }
1302 return gsl_linalg_COD_lssolve( QRZT.
get(), tau_Q.
get(), tau_Z.
get(), perm.
get(),
1323 return gsl_linalg_COD_lssolve2( lambda, QRZT.
get(), tau_Q.
get(), tau_Z.
get(), perm.
get(),
1339 return gsl_linalg_COD_unpack( QRZT.
get(), tau_Q.
get(), tau_Z.
get(), rank, Q.
get(),
1352 return gsl_linalg_COD_matZ( QRZT.
get(), tau_Z.
get(), rank, A.
get(), work.
get() ); }
1373 return gsl_linalg_LQ_QTvec( LQ.
get(), tau.
get(), v.
get() ); }
1388 return gsl_linalg_cholesky_solve_mat( cholesky.
get(), B.
get(), X.
get() ); }
1396 return gsl_linalg_cholesky_svx_mat( cholesky.
get(), X.
get() ); }
1404 return gsl_linalg_cholesky_scale( A.
get(), S.
get() ); }
1412 return gsl_linalg_cholesky_scale_apply( A.
get(), S.
get() ); }
1420 return gsl_linalg_cholesky_decomp2( A.
get(), S.
get() ); }
1429 return gsl_linalg_cholesky_svx2( LLT.
get(), S.
get(), x.
get() ); }
1440 return gsl_linalg_cholesky_solve2( LLT.
get(), S.
get(),
b.get(), x.
get() ); }
1449 return gsl_linalg_cholesky_rcond( LLT.
get(), &
rcond, work.
get() ); }
1457 return gsl_linalg_pcholesky_decomp( A.
get(), p.
get() ); }
1468 return gsl_linalg_pcholesky_solve( LDLT.
get(), p.
get(),
b.get(), x.
get() ); }
1477 return gsl_linalg_pcholesky_svx( LDLT.
get(), p.
get(), x.
get() ); }
1486 return gsl_linalg_pcholesky_decomp2( A.
get(), p.
get(), S.
get() ); }
1498 return gsl_linalg_pcholesky_solve2( LDLT.
get(), p.
get(), S.
get(),
b.get(), x.
get() ); }
1509 return gsl_linalg_pcholesky_svx2( LDLT.
get(), p.
get(), S.
get(), x.
get() ); }
1518 return gsl_linalg_pcholesky_invert( LDLT.
get(), p.
get(), Ainv.
get() ); }
1529 return gsl_linalg_pcholesky_rcond( LDLT.
get(), p.
get(), &
rcond, work.
get() ); }
1538 return gsl_linalg_mcholesky_decomp( A.
get(), p.
get(),
E.get() ); }
1549 return gsl_linalg_mcholesky_solve( LDLT.
get(), p.
get(),
b.get(), x.
get() ); }
1558 return gsl_linalg_mcholesky_svx( LDLT.
get(), p.
get(), x.
get() ); }
1569 return gsl_linalg_mcholesky_rcond( LDLT.
get(), p.
get(), &
rcond, work.
get() ); }
1578 return gsl_linalg_mcholesky_invert( LDLT.
get(), p.
get(), Ainv.
get() ); }
1587 return gsl_linalg_cholesky_band_decomp( A.
get() ); }
1596 return gsl_linalg_cholesky_band_solve( LLT.
get(),
b.get(), x.
get() ); }
1604 return gsl_linalg_cholesky_band_svx( LLT.
get(), x.
get() ); }
1612 return gsl_linalg_cholesky_band_invert( LLT.
get(), Ainv.
get() ); }
1620 return gsl_linalg_cholesky_band_unpack( LLT.
get(), L.
get() ); }
1629 return gsl_linalg_cholesky_band_rcond( LLT.
get(), &
rcond, work.
get() ); }
1644 return gsl_linalg_ldlt_solve( LDLT.
get(),
b.get(), x.
get() ); }
1652 return gsl_linalg_ldlt_svx( LDLT.
get(), x.
get() ); }
1661 return gsl_linalg_ldlt_rcond( LDLT.
get(), &
rcond, work.
get() ); }
1669 return gsl_linalg_ldlt_band_decomp( A.
get() ); }
1678 return gsl_linalg_ldlt_band_solve( LDLT.
get(),
b.get(), x.
get() ); }
1686 return gsl_linalg_ldlt_band_svx( LDLT.
get(), x.
get() ); }
1695 return gsl_linalg_ldlt_band_unpack( LDLT.
get(), L.
get(),
D.get() ); }
1704 return gsl_linalg_ldlt_band_rcond( LDLT.
get(),
rcond, work.
get() ); }
1713 return gsl_linalg_ldlt_band_rcond( LDLT.
get(), &
rcond, work.
get() ); }
1723 return gsl_linalg_tri_rcond( Uplo, A.
get(), &
rcond, work.
get() ); }
1732 return gsl_linalg_tri_upper_rcond( A.
get(), &
rcond, work.
get() ); }
1741 return gsl_linalg_tri_lower_rcond( A.
get(), &
rcond, work.
get() ); }
1752 inline int invnorm1(
size_t const N,
1753 int (*Ainvx)(CBLAS_TRANSPOSE_t TransA, gsl_vector* x,
void* params),
1754 void* params,
double* Ainvnorm,
vector& work ){
1755 return gsl_linalg_invnorm1( N, Ainvx, params, Ainvnorm, work.
get() ); }
1776 return gsl_linalg_tri_upper_unit_invert( T.
get() ); }
1783 return gsl_linalg_tri_lower_unit_invert( T.
get() ); }
1792 return gsl_linalg_tri_invert( Uplo, Diag, T.
get() ); }
1801 return gsl_linalg_complex_tri_invert( Uplo, Diag, T.
get() ); }
1820 return gsl_linalg_complex_tri_LHL( L.
get() ); }
1827 return gsl_linalg_complex_tri_UL( LU.
get() ); }
1835 inline void givens(
double const a,
double const b,
double& c,
double& s ){
1836 return gsl_linalg_givens(
a,
b, &c, &s ); }
1847 return gsl_linalg_givens_gv( v.
get(), i, j, c, s ); }
1859 return gsl_linalg_LU_band_decomp( M, lb, ub, AB.
get(), piv.
get() ); }
1873 return gsl_linalg_LU_band_solve( lb, ub, LUB.
get(), piv.
get(),
b.get(), x.
get() ); }
1885 return gsl_linalg_LU_band_svx( lb, ub, LUB.
get(), piv.
get(), x.
get() ); }
1900 return gsl_linalg_LU_band_unpack( M, lb, ub, LUB.
get(), piv.
get(), L.
get(), U.
get() ); }
1908 return gsl_linalg_QR_decomp_old( A.
get(), tau.
get() ); }
1916 return gsl_linalg_complex_QR_decomp( A.
get(), tau.
get() ); }
1924 return gsl_linalg_complex_QR_decomp_r( A.
get(), T.
get() ); }
1935 return gsl_linalg_complex_QR_solve( QR.
get(), tau.
get(),
b.get(), x.
get() ); }
1946 return gsl_linalg_complex_QR_solve_r( QR.
get(), T.
get(),
b.get(), x.
get() ); }
1956 return gsl_linalg_complex_QR_svx( QR.
get(), tau.
get(), x.
get() ); }
1969 return gsl_linalg_complex_QR_lssolve( QR.
get(), tau.
get(),
b.get(), x.
get(),
1983 return gsl_linalg_complex_QR_lssolve_r( QR.
get(), T.
get(),
b.get(), x.
get(),
1994 return gsl_linalg_complex_QR_QHvec( QR.
get(), tau.
get(), v.
get() ); }
2005 return gsl_linalg_complex_QR_QHvec_r( QR.
get(), T.
get(),
b.get(), work.
get() ); }
2015 return gsl_linalg_complex_QR_Qvec( QR.
get(), tau.
get(), v.
get() ); }
2026 return gsl_linalg_complex_QR_unpack( QR.
get(), tau.
get(), Q.
get(), R.
get() ); }
2037 return gsl_linalg_complex_QR_unpack_r( QR.
get(), T.
get(), Q.
get(), R.
get() ); }
2049 return gsl_linalg_QR_band_decomp_L2( M, p, q, AB.
get(), tau.
get() ); }
2062 return gsl_linalg_QR_band_unpack_L2( p, q, QRB.
get(), tau.
get(), Q.
get(), R.
get() ); }
2072 return gsl_linalg_QR_UD_decomp( U.
get(),
D.get(), Y.
get(), T.
get() ); }
2085 return gsl_linalg_QR_UD_lssolve( R.
get(), Y.
get(), T.
get(),
b.get(), x.
get(),
2095 return gsl_linalg_QR_UR_decomp( S.
get(), A.
get(), T.
get() ); }
2104 return gsl_linalg_QR_UU_decomp( U.
get(), S.
get(), T.
get() ); }
2117 return gsl_linalg_QR_UU_lssolve( R.
get(), Y.
get(), T.
get(),
b.get(), x.
get(),
2128 return gsl_linalg_QR_UU_QTvec( Y.
get(), T.
get(),
b.get(), work.
get() ); }
2137 return gsl_linalg_QR_UZ_decomp( S.
get(), A.
get(), T.
get() ); }
2145 return gsl_linalg_QL_decomp( A.
get(), tau.
get() ); }
2155 return gsl_linalg_QL_unpack( QL.
get(), tau.
get(), Q.
get(), L.
get() ); }
This class handles complex numbers.
gsl_complex & get()
Get the base class object.
This class handles matrix_complex objects as shared handles.
gsl_matrix_complex * get()
Get the gsl_matrix_complex.
This class handles matrix objects as shared handles.
gsl_matrix * get()
Get the gsl_matrix.
This class handles GSL permutation objects.
gsl_permutation * get()
Get the gsl_permutation.
This class handles vector_complex objects as shared handles.
gsl_vector_complex * get()
Get the gsl_vector_complex.
This class handles vector_uint objects as shared handles.
gsl_vector_uint * get()
Get the gsl_vector_uint.
This class handles vector objects as shared handles.
gsl_vector * get()
Get the gsl_vector.
complex inverse(complex const &a)
C++ version of gsl_complex_inverse().
double LU_lndet(matrix &LU)
C++ version of gsl_linalg_LU_lndet().
int QRPT_lssolve(matrix const &QR, vector const &tau, permutation const &p, vector const &b, vector &x, vector &residual)
C++ version of gsl_linalg_QRPT_lssolve().
int cholesky_svx2(matrix const &LLT, vector const &S, vector &x)
C++ version of gsl_linalg_cholesky_svx2().
int mcholesky_rcond(matrix const &LDLT, permutation const &p, double &rcond, vector &work)
C++ version of gsl_linalg_mcholesky_rcond().
int LU_solve(matrix const &LU, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_LU_solve().
int complex_QR_unpack(matrix_complex const &QR, vector_complex const &tau, matrix_complex &Q, matrix_complex &R)
C++ version of gsl_linalg_complex_QR_unpack().
int pcholesky_decomp(matrix &A, permutation &p)
C++ version of gsl_linalg_pcholesky_decomp().
int cholesky_scale_apply(matrix &A, vector const &S)
C++ version of gsl_linalg_cholesky_scale_apply().
int balance_matrix(matrix &A, vector &D)
C++ version of gsl_linalg_balance_matrix().
int LU_band_svx(size_t const lb, size_t const ub, matrix const &LUB, vector_uint const &piv, vector &x)
C++ version of gsl_linalg_LU_band_svx().
int bidiag_decomp(matrix &A, vector &tau_U, vector &tau_V)
C++ version of gsl_linalg_bidiag_decomp().
int QR_Qvec(matrix const &QR, vector const &tau, vector &v)
C++ version of gsl_linalg_QR_Qvec().
int cholesky_band_invert(matrix const &LLT, matrix &Ainv)
C++ version of gsl_linalg_cholesky_band_invert().
int complex_QR_decomp(matrix_complex &A, vector_complex &tau)
C++ version of gsl_linalg_complex_QR_decomp().
int complex_QR_QHvec_r(matrix_complex const &QR, matrix_complex const &T, vector_complex &b, vector_complex &work)
C++ version of gsl_linalg_complex_QR_QHvec_r().
int complex_QR_unpack_r(matrix_complex const &QR, matrix_complex const &T, matrix_complex &Q, matrix_complex &R)
C++ version of gsl_linalg_complex_QR_unpack_r().
int cholesky_solve(matrix const &cholesky, vector const &b, vector &x)
C++ version of gsl_linalg_cholesky_solve().
int ldlt_svx(matrix const &LDLT, vector &x)
C++ version of gsl_linalg_ldlt_svx().
int LU_decomp(matrix &A, permutation &p, int &signum)
C++ version of gsl_linalg_LU_decomp().
int balance_columns(matrix &A, vector &D)
C++ version of gsl_linalg_balance_columns().
int R_svx(matrix const &R, vector &x)
C++ version of gsl_linalg_R_svx().
int pcholesky_solve(matrix const &LDLT, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_pcholesky_solve().
int QRPT_svx(matrix const &QR, vector const &tau, permutation const &p, vector &x)
C++ version of gsl_linalg_QRPT_svx().
int solve_symm_tridiag(vector const &diag, vector const &offdiag, vector const &b, vector &x)
C++ version of gsl_linalg_solve_symm_tridiag().
complex complex_LU_det(matrix_complex &LU, int signum)
C++ version of gsl_linalg_complex_LU_det().
int cholesky_svx_mat(matrix const &cholesky, matrix &X)
C++ version of gsl_linalg_cholesky_svx_mat().
int cholesky_band_unpack(matrix const &LLT, matrix &L)
C++ version of gsl_linalg_cholesky_band_unpack().
int COD_lssolve(matrix const &QRZT, vector const &tau_Q, vector const &tau_Z, permutation const &perm, size_t const rank, vector const &b, vector &x, vector &residual)
C++ version of gsl_linalg_COD_lssolve().
complex complex_LU_sgndet(matrix_complex &LU, int signum)
C++ version of gsl_linalg_complex_LU_sgndet().
int ldlt_band_rcond(matrix const &LDLT, double *rcond, vector &work)
C++ version of gsl_linalg_ldlt_band_rcond().
int hessenberg_set_zero(matrix &H)
C++ version of gsl_linalg_hessenberg_set_zero().
int mcholesky_svx(matrix const &LDLT, permutation const &p, vector &x)
C++ version of gsl_linalg_mcholesky_svx().
int tri_invert(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix &T)
C++ version of gsl_linalg_tri_invert().
int complex_QR_QHvec(matrix_complex const &QR, vector_complex const &tau, vector_complex &v)
C++ version of gsl_linalg_complex_QR_QHvec().
int COD_decomp_e(matrix &A, vector &tau_Q, vector &tau_Z, permutation &p, double tol, size_t &rank, vector &work)
C++ version of gsl_linalg_COD_decomp_e().
int complex_QR_lssolve_r(matrix_complex const &QR, matrix_complex const &T, vector_complex const &b, vector_complex &x, vector_complex &work)
C++ version of gsl_linalg_complex_QR_lssolve_r().
int QR_Rsolve(matrix const &QR, vector const &b, vector &x)
C++ version of gsl_linalg_QR_Rsolve().
int LQ_decomp(matrix &A, vector &tau)
C++ version of gsl_linalg_LQ_decomp().
int bidiag_unpack_B(matrix const &A, vector &diag, vector &superdiag)
C++ version of gsl_linalg_bidiag_unpack_B().
void givens_gv(vector &v, size_t const i, size_t const j, double const c, double const s)
C++ version of gsl_linalg_givens_gv().
int QRPT_rcond(matrix const &QR, double &rcond, vector &work)
C++ version of gsl_linalg_QRPT_rcond().
int QR_UU_decomp(matrix &U, matrix &S, matrix &T)
C++ version of gsl_linalg_QR_UU_decomp().
int PTLQ_decomp2(matrix const &A, matrix &q, matrix &r, vector &tau, permutation &p, int *signum, vector &norm)
C++ version of gsl_linalg_PTLQ_decomp2().
int QL_unpack(matrix const &QL, vector const &tau, matrix &Q, matrix &L)
C++ version of gsl_linalg_QL_unpack().
int SV_decomp(matrix &A, matrix &V, vector &S, vector &work)
C++ version of gsl_linalg_SV_decomp().
int householder_hm1(double tau, matrix &A)
C++ version of gsl_linalg_householder_hm1().
int cholesky_band_solve(matrix const &LLT, vector const &b, vector &x)
C++ version of gsl_linalg_cholesky_band_solve().
int cholesky_solve_mat(matrix const &cholesky, matrix const &B, matrix &X)
C++ version of gsl_linalg_cholesky_solve_mat().
int LQ_svx_T(matrix const &LQ, vector const &tau, vector &x)
C++ version of gsl_linalg_LQ_svx_T().
int QR_unpack(matrix const &QR, vector const &tau, matrix &Q, matrix &R)
C++ version of gsl_linalg_QR_unpack().
int complex_cholesky_solve(matrix_complex const &cholesky, vector_complex const &b, vector_complex &x)
C++ version of gsl_linalg_complex_cholesky_solve().
int QRPT_decomp2(matrix const &A, matrix &q, matrix &r, vector &tau, permutation &p, int &signum, vector &norm)
C++ version of gsl_linalg_QRPT_decomp2().
int LU_band_unpack(size_t const M, size_t const lb, size_t const ub, matrix const &LUB, vector_uint const &piv, matrix &L, matrix &U)
C++ version of gsl_linalg_LU_band_unpack().
int cholesky_decomp_unit(matrix &A, vector &D)
C++ version of gsl_linalg_cholesky_decomp_unit().
int LQ_vecQ(matrix const &LQ, vector const &tau, vector &v)
C++ version of gsl_linalg_LQ_vecQ().
int QRPT_update(matrix &Q, matrix &R, permutation const &p, vector &u, vector const &v)
C++ version of gsl_linalg_QRPT_update().
int complex_tri_invert(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix_complex &T)
C++ version of gsl_linalg_complex_tri_invert().
int solve_cyc_tridiag(vector const &diag, vector const &abovediag, vector const &belowdiag, vector const &b, vector &x)
C++ version of gsl_linalg_solve_cyc_tridiag().
int complex_LU_invx(matrix_complex &LU, permutation const &p)
C++ version of gsl_linalg_complex_LU_invx().
int LQ_solve_T(matrix const &LQ, vector const &tau, vector const &b, vector &x)
C++ version of gsl_linalg_LQ_solve_T().
int complex_LU_solve(matrix_complex const &LU, permutation const &p, vector_complex const &b, vector_complex &x)
C++ version of gsl_linalg_complex_LU_solve().
int cholesky_decomp1(matrix &A)
C++ version of gsl_linalg_cholesky_decomp1().
int QR_UR_decomp(matrix &S, matrix &A, matrix &T)
C++ version of gsl_linalg_QR_UR_decomp().
int tri_UL(matrix &LU)
C++ version of gsl_linalg_tri_UL().
int pcholesky_rcond(matrix const &LDLT, permutation const &p, double &rcond, vector &work)
C++ version of gsl_linalg_pcholesky_rcond().
int HH_svx(matrix &A, vector &x)
C++ version of gsl_linalg_HH_svx().
int complex_LU_refine(matrix_complex const &A, matrix_complex const &LU, permutation const &p, vector_complex const &b, vector_complex &x, vector_complex &residual)
C++ version of gsl_linalg_complex_LU_refine().
int QR_unpack_r(matrix const &QR, matrix const &T, matrix &Q, matrix &R)
C++ version of gsl_linalg_QR_unpack_r().
int ldlt_band_unpack(matrix const &LDLT, matrix &L, vector &D)
C++ version of gsl_linalg_ldlt_band_unpack().
int tri_upper_invert(matrix &T)
C++ version of gsl_linalg_tri_upper_invert().
int LQ_Lsvx_T(matrix const &LQ, vector &x)
C++ version of gsl_linalg_LQ_Lsvx_T().
int cholesky_decomp2(matrix &A, vector &S)
C++ version of gsl_linalg_cholesky_decomp2().
int SV_decomp_mod(matrix &A, matrix &X, matrix &V, vector &S, vector &work)
C++ version of gsl_linalg_SV_decomp_mod().
int QR_UU_lssolve(matrix const &R, matrix const &Y, matrix const &T, vector const &b, vector &x, vector &work)
C++ version of gsl_linalg_QR_UU_lssolve().
int cholesky_svx(matrix const &cholesky, vector &x)
C++ version of gsl_linalg_cholesky_svx().
int COD_matZ(matrix const &QRZT, vector const &tau_Z, size_t const rank, matrix &A, vector &work)
C++ version of gsl_linalg_COD_matZ().
int QRPT_lssolve2(matrix const &QR, vector const &tau, permutation const &p, vector const &b, size_t const rank, vector &x, vector &residual)
C++ version of gsl_linalg_QRPT_lssolve2().
int QR_decomp(matrix &A, vector &tau)
C++ version of gsl_linalg_QR_decomp().
int LU_band_solve(size_t const lb, size_t const ub, matrix const &LUB, vector_uint const &piv, vector const &b, vector &x)
C++ version of gsl_linalg_LU_band_solve().
int LQ_Lsolve_T(matrix const &LQ, vector const &b, vector &x)
C++ version of gsl_linalg_LQ_Lsolve_T().
int QR_QTmat(matrix const &QR, vector const &tau, matrix &A)
C++ version of gsl_linalg_QR_QTmat().
int pcholesky_decomp2(matrix &A, permutation &p, vector &S)
C++ version of gsl_linalg_pcholesky_decomp2().
int bidiag_unpack(matrix const &A, vector const &tau_U, matrix &U, vector const &tau_V, matrix &V, vector &diag, vector &superdiag)
C++ version of gsl_linalg_bidiag_unpack().
int COD_decomp(matrix &A, vector &tau_Q, vector &tau_Z, permutation &p, size_t &rank, vector &work)
C++ version of gsl_linalg_COD_decomp().
int complex_QR_Qvec(matrix_complex const &QR, vector_complex const &tau, vector_complex &v)
C++ version of gsl_linalg_complex_QR_Qvec().
int COD_unpack(matrix const &QRZT, vector const &tau_Q, vector const &tau_Z, size_t const rank, matrix &Q, matrix &R, matrix &Z)
C++ version of gsl_linalg_COD_unpack().
double householder_transform2(double &alpha, vector &v)
C++ version of gsl_linalg_householder_transform2().
int SV_decomp_jacobi(matrix &A, matrix &Q, vector &S)
C++ version of gsl_linalg_SV_decomp_jacobi().
int complex_cholesky_svx(matrix_complex const &cholesky, vector_complex &x)
C++ version of gsl_linalg_complex_cholesky_svx().
int tri_upper_unit_invert(matrix &T)
C++ version of gsl_linalg_tri_upper_unit_invert().
int tri_LTL(matrix &L)
C++ version of gsl_linalg_tri_LTL().
int symmtd_unpack(matrix const &A, vector const &tau, matrix &Q, vector &diag, vector &subdiag)
C++ version of gsl_linalg_symmtd_unpack().
int hermtd_unpack_T(matrix_complex const &A, vector &diag, vector &subdiag)
C++ version of gsl_linalg_hermtd_unpack_T().
int QR_QTvec(matrix const &QR, vector const &tau, vector &v)
C++ version of gsl_linalg_QR_QTvec().
int QR_band_unpack_L2(size_t const p, size_t const q, matrix const &QRB, vector const &tau, matrix &Q, matrix &R)
C++ version of gsl_linalg_QR_band_unpack_L2().
int hessenberg_unpack_accum(matrix &H, vector &tau, matrix &U)
C++ version of gsl_linalg_hessenberg_unpack_accum().
int hermtd_decomp(matrix_complex &A, vector_complex &tau)
C++ version of gsl_linalg_hermtd_decomp().
int QL_decomp(matrix &A, vector &tau)
C++ version of gsl_linalg_QL_decomp().
int PTLQ_svx_T(matrix const &LQ, vector const &tau, permutation const &p, vector &x)
C++ version of gsl_linalg_PTLQ_svx_T().
complex complex_householder_transform(vector_complex &v)
C++ version of gsl_linalg_complex_householder_transform().
int QR_UU_QTvec(matrix const &Y, matrix const &T, vector &b, vector &work)
C++ version of gsl_linalg_QR_UU_QTvec().
int cholesky_band_rcond(matrix const &LLT, double &rcond, vector &work)
C++ version of gsl_linalg_cholesky_band_rcond().
int cholesky_decomp(matrix &A)
C++ version of gsl_linalg_cholesky_decomp().
int QR_update(matrix &Q, matrix &R, vector &w, vector const &v)
C++ version of gsl_linalg_QR_update().
int LU_refine(matrix const &A, matrix const &LU, permutation const &p, vector const &b, vector &x, vector &residual)
C++ version of gsl_linalg_LU_refine().
int householder_hv(double tau, vector const &v, vector &w)
C++ version of gsl_linalg_householder_hv().
int LQ_LQsolve(matrix &Q, matrix &L, vector const &b, vector &x)
C++ version of gsl_linalg_LQ_LQsolve().
double LU_det(matrix &LU, int signum)
C++ version of gsl_linalg_LU_det().
int complex_QR_solve(matrix_complex const &QR, vector_complex const &tau, vector_complex const &b, vector_complex &x)
C++ version of gsl_linalg_complex_QR_solve().
int complex_householder_mh(complex &tau, vector_complex const &v, matrix_complex &A)
C++ version of gsl_linalg_complex_householder_mh().
int ldlt_rcond(matrix const &LDLT, double &rcond, vector &work)
C++ version of gsl_linalg_ldlt_rcond().
int complex_tri_LHL(matrix_complex &L)
C++ version of gsl_linalg_complex_tri_LHL().
int complex_LU_svx(matrix_complex const &LU, permutation const &p, vector_complex &x)
C++ version of gsl_linalg_complex_LU_svx().
int QR_rcond(matrix const &QR, double &rcond, vector &work)
C++ version of gsl_linalg_QR_rcond().
int complex_QR_solve_r(matrix_complex const &QR, matrix_complex const &T, vector_complex const &b, vector_complex &x)
C++ version of gsl_linalg_complex_QR_solve_r().
int tri_rcond(CBLAS_UPLO_t Uplo, matrix const &A, double &rcond, vector &work)
C++ version of gsl_linalg_tri_rcond().
int mcholesky_decomp(matrix &A, permutation &p, vector &E)
C++ version of gsl_linalg_mcholesky_decomp().
int pcholesky_invert(matrix const &LDLT, permutation const &p, matrix &Ainv)
C++ version of gsl_linalg_pcholesky_invert().
int LQ_update(matrix &Q, matrix &R, vector const &v, vector &w)
C++ version of gsl_linalg_LQ_update().
int symmtd_unpack_T(matrix const &A, vector &diag, vector &subdiag)
C++ version of gsl_linalg_symmtd_unpack_T().
int PTLQ_LQsolve_T(matrix const &Q, matrix const &L, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_PTLQ_LQsolve_T().
int ldlt_band_solve(matrix const &LDLT, vector const &b, vector &x)
C++ version of gsl_linalg_ldlt_band_solve().
int QRPT_Rsolve(matrix const &QR, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_QRPT_Rsolve().
int mcholesky_invert(matrix const &LDLT, permutation const &p, matrix &Ainv)
C++ version of gsl_linalg_mcholesky_invert().
int SV_solve(matrix const &U, matrix const &Q, vector const &S, vector const &b, vector &x)
C++ version of gsl_linalg_SV_solve().
int QRPT_QRsolve(matrix const &Q, matrix const &R, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_QRPT_QRsolve().
int symmtd_decomp(matrix &A, vector &tau)
C++ version of gsl_linalg_symmtd_decomp().
int LU_invx(matrix &LU, permutation const &p)
C++ version of gsl_linalg_LU_invx().
int PTLQ_solve_T(matrix const &QR, vector const &tau, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_PTLQ_solve_T().
int QR_Rsvx(matrix const &QR, vector &x)
C++ version of gsl_linalg_QR_Rsvx().
int QRPT_decomp(matrix &A, vector &tau, permutation &p, int &signum, vector &norm)
C++ version of gsl_linalg_QRPT_decomp().
int cholesky_band_svx(matrix const &LLT, vector &x)
C++ version of gsl_linalg_cholesky_band_svx().
int complex_householder_hv(complex &tau, vector_complex const &v, vector_complex &w)
C++ version of gsl_linalg_complex_householder_hv().
int R_solve(matrix const &R, vector const &b, vector &x)
C++ version of gsl_linalg_R_solve().
int complex_LU_invert(matrix_complex const &LU, permutation const &p, matrix_complex &inverse)
C++ version of gsl_linalg_complex_LU_invert().
int hesstri_decomp(matrix &A, matrix &B, matrix &U, matrix &V, vector &work)
C++ version of gsl_linalg_hesstri_decomp().
int QR_UD_decomp(matrix &U, vector const &D, matrix &Y, matrix &T)
C++ version of gsl_linalg_QR_UD_decomp().
int QR_matQ(matrix const &QR, vector const &tau, matrix &A)
C++ version of gsl_linalg_QR_matQ().
size_t QRPT_rank(matrix const &QR, double const tol)
C++ version of gsl_linalg_QRPT_rank().
int QR_solve(matrix const &QR, vector const &tau, vector const &b, vector &x)
C++ version of gsl_linalg_QR_solve().
int QRPT_solve(matrix const &QR, vector const &tau, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_QRPT_solve().
int QR_decomp_old(matrix &A, vector &tau)
C++ version of gsl_linalg_QR_decomp_old().
int cholesky_invert(matrix &cholesky)
C++ version of gsl_linalg_cholesky_invert().
int HH_solve(matrix &A, vector const &b, vector &x)
C++ version of gsl_linalg_HH_solve().
int L_solve_T(matrix const &L, vector const &b, vector &x)
C++ version of gsl_linalg_L_solve_T().
int complex_householder_hm(complex &tau, vector_complex const &v, matrix_complex &A)
C++ version of gsl_linalg_complex_householder_hm().
int hessenberg_unpack(matrix &H, vector &tau, matrix &U)
C++ version of gsl_linalg_hessenberg_unpack().
int solve_symm_cyc_tridiag(vector const &diag, vector const &offdiag, vector const &b, vector &x)
C++ version of gsl_linalg_solve_symm_cyc_tridiag().
int pcholesky_svx(matrix const &LDLT, permutation const &p, vector &x)
C++ version of gsl_linalg_pcholesky_svx().
int mcholesky_solve(matrix const &LDLT, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_mcholesky_solve().
int cholesky_rcond(matrix const &LLT, double &rcond, vector &work)
C++ version of gsl_linalg_cholesky_rcond().
int SV_leverage(matrix const &U, vector &h)
C++ version of gsl_linalg_SV_leverage().
int complex_QR_lssolve(matrix_complex const &QR, vector_complex const &tau, vector_complex const &b, vector_complex &x, vector_complex &residual)
C++ version of gsl_linalg_complex_QR_lssolve().
int LQ_unpack(matrix const &LQ, vector const &tau, matrix &Q, matrix &L)
C++ version of gsl_linalg_LQ_unpack().
int QR_UZ_decomp(matrix &S, matrix &A, matrix &T)
C++ version of gsl_linalg_QR_UZ_decomp().
int PTLQ_decomp(matrix &A, vector &tau, permutation &p, int &signum, vector &norm)
C++ version of gsl_linalg_PTLQ_decomp().
int tri_lower_rcond(matrix const &A, double &rcond, vector &work)
C++ version of gsl_linalg_tri_lower_rcond().
int LQ_vecQT(matrix const &LQ, vector const &tau, vector &v)
C++ version of gsl_linalg_LQ_vecQT().
int complex_QR_svx(matrix_complex const &QR, vector_complex const &tau, vector_complex &x)
C++ version of gsl_linalg_complex_QR_svx().
int PTLQ_Lsolve_T(matrix const &LQ, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_PTLQ_Lsolve_T().
int LU_svx(matrix const &LU, permutation const &p, vector &x)
C++ version of gsl_linalg_LU_svx().
int PTLQ_update(matrix &Q, matrix &L, permutation const &p, vector const &v, vector &w)
C++ version of gsl_linalg_PTLQ_update().
int QR_lssolve(matrix const &QR, vector const &tau, vector const &b, vector &x, vector &residual)
C++ version of gsl_linalg_QR_lssolve().
double householder_transform(vector &v)
C++ version of gsl_linalg_householder_transform().
int ldlt_decomp(matrix &A)
C++ version of gsl_linalg_ldlt_decomp().
int solve_tridiag(vector const &diag, vector const &abovediag, vector const &belowdiag, vector const &b, vector &x)
C++ version of gsl_linalg_solve_tridiag().
int PTLQ_Lsvx_T(matrix const &LQ, permutation const &p, vector &x)
C++ version of gsl_linalg_PTLQ_Lsvx_T().
int hessenberg_decomp(matrix &A, vector &tau)
C++ version of gsl_linalg_hessenberg_decomp().
int complex_cholesky_decomp(matrix_complex &A)
C++ version of gsl_linalg_complex_cholesky_decomp().
int QR_UD_lssolve(matrix const &R, matrix const &Y, matrix const &T, vector const &b, vector &x, vector &work)
C++ version of gsl_linalg_QR_UD_lssolve().
int complex_tri_UL(matrix_complex &LU)
C++ version of gsl_linalg_complex_tri_UL().
int pcholesky_svx2(matrix const &LDLT, permutation const &p, vector const &S, vector &x)
C++ version of gsl_linalg_pcholesky_svx2().
int LU_invert(matrix const &LU, permutation const &p, matrix &inverse)
C++ version of gsl_linalg_LU_invert().
int bidiag_unpack2(matrix &A, vector &tau_U, vector &tau_V, matrix &V)
C++ version of gsl_linalg_bidiag_unpack2().
int householder_mh(double tau, vector const &v, matrix &A)
C++ version of gsl_linalg_householder_mh().
int householder_left(double const tau, vector const &v, matrix &A, vector &work)
C++ version of gsl_linalg_householder_left().
double complex_LU_lndet(matrix_complex &LU)
C++ version of gsl_linalg_complex_LU_lndet().
int LU_band_decomp(size_t const M, size_t const lb, size_t const ub, matrix &AB, vector_uint &piv)
C++ version of gsl_linalg_LU_band_decomp().
int householder_right(double const tau, vector const &v, matrix &A, vector &work)
C++ version of gsl_linalg_householder_right().
int QR_QTmat_r(matrix const &QR, matrix const &T, matrix &B, matrix &work)
C++ version of gsl_linalg_QR_QTmat_r().
int QR_band_decomp_L2(size_t const M, size_t const p, size_t const q, matrix &AB, vector &tau)
C++ version of gsl_linalg_QR_band_decomp_L2().
int hessenberg_submatrix(matrix &M, matrix &A, size_t top, vector &tau)
C++ version of gsl_linalg_hessenberg_submatrix().
int complex_LU_decomp(matrix_complex &A, permutation &p, int &signum)
C++ version of gsl_linalg_complex_LU_decomp().
int QR_svx(matrix const &QR, vector const &tau, vector &x)
C++ version of gsl_linalg_QR_svx().
int householder_hm(double tau, vector const &v, matrix &A)
C++ version of gsl_linalg_householder_hm().
int cholesky_solve2(matrix const &LLT, vector const &S, vector const &b, vector &x)
C++ version of gsl_linalg_cholesky_solve2().
int QR_QRsolve(matrix &Q, matrix &R, vector const &b, vector &x)
C++ version of gsl_linalg_QR_QRsolve().
int cholesky_band_decomp(matrix &A)
C++ version of gsl_linalg_cholesky_band_decomp().
int complex_QR_decomp_r(matrix_complex &A, matrix_complex &T)
C++ version of gsl_linalg_complex_QR_decomp_r().
int ldlt_solve(matrix const &LDLT, vector const &b, vector &x)
C++ version of gsl_linalg_ldlt_solve().
int LQ_lssolve(matrix const &LQ, vector const &tau, vector const &b, vector &x, vector &residual)
C++ version of gsl_linalg_LQ_lssolve().
void givens(double const a, double const b, double &c, double &s)
C++ version of gsl_linalg_givens().
int QRPT_Rsvx(matrix const &QR, permutation const &p, vector &x)
C++ version of gsl_linalg_QRPT_Rsvx().
int LU_sgndet(matrix &lu, int signum)
C++ version of gsl_linalg_LU_sgndet().
int COD_lssolve2(double const lambda, matrix const &QRZT, vector const &tau_Q, vector const &tau_Z, permutation const &perm, size_t const rank, vector const &b, vector &x, vector &residual, matrix &S, vector &work)
C++ version of gsl_linalg_COD_lssolve2().
int ldlt_band_svx(matrix const &LDLT, vector &x)
C++ version of gsl_linalg_ldlt_band_svx().
int balance_accum(matrix &A, vector &D)
C++ version of gsl_linalg_balance_accum().
int cholesky_scale(matrix const &A, vector &S)
C++ version of gsl_linalg_cholesky_scale().
int tri_upper_rcond(matrix const &A, double &rcond, vector &work)
C++ version of gsl_linalg_tri_upper_rcond().
int LQ_lssolve_T(matrix const &LQ, vector const &tau, vector const &b, vector &x, vector &residual)
C++ version of gsl_linalg_LQ_lssolve_T().
int tri_lower_unit_invert(matrix &T)
C++ version of gsl_linalg_tri_lower_unit_invert().
int hermtd_unpack(matrix_complex const &A, vector_complex const &tau, matrix_complex &U, vector &diag, vector &sudiag)
C++ version of gsl_linalg_hermtd_unpack().
int LQ_QTvec(matrix const &LQ, vector const &tau, vector &v)
C++ version of gsl_linalg_LQ_QTvec().
int pcholesky_solve2(matrix const &LDLT, permutation const &p, vector const &S, vector const &b, vector &x)
C++ version of gsl_linalg_pcholesky_solve2().
int complex_householder_left(complex const tau, vector_complex const &v, matrix_complex &A, vector_complex &work)
C++ version of gsl_linalg_complex_householder_left().
int tri_lower_invert(matrix &T)
C++ version of gsl_linalg_tri_lower_invert().
int ldlt_band_decomp(matrix &A)
C++ version of gsl_linalg_ldlt_band_decomp().
int rcond(double &rcond, workspace const &w)
C++ version of gsl_multifit_nlinear_rcond().
vector residual(workspace const &w)
C++ version of gsl_multifit_nlinear_residual().
double D(double phi, double k, mode_t mode)
C++ version of gsl_sf_ellint_D().
double b(int order, double qq)
C++ version of gsl_sf_mathieu_b().
double a(int order, double qq)
C++ version of gsl_sf_mathieu_a().
The gsl package creates an interface to the GNU Scientific Library for C++.