ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
linalg.hpp
Go to the documentation of this file.
1/*
2 * $Id: linalg.hpp 308 2013-11-17 16:25:51Z jdl3 $
3 * Copyright (C) 2010, 2018, 2019, 2020, 2021 John D Lamb
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or (at
8 * your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful, but
11 * WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 * General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
18 */
19
20#ifndef CCGSL_LINALG_HPP
21#define CCGSL_LINALG_HPP
22
23#include<cmath>
24#include<gsl/gsl_linalg.h>
25#include"vector.hpp"
26#include"matrix.hpp"
27#include"complex.hpp"
28#include"permutation.hpp"
29#include"vector_complex.hpp"
30#include"vector_uint.hpp"
31#include"matrix_complex.hpp"
32
33namespace gsl {
37 namespace linalg {
43 inline double householder_transform( vector& v ){
44 return gsl_linalg_householder_transform( v.get() ); }
51 inline double householder_transform2( double& alpha, vector& v ){
52 return gsl_linalg_householder_transform2( &alpha, v.get() ); }
59 return gsl_linalg_complex_householder_transform( v.get() ); }
67 inline int householder_hm( double tau, vector const& v, matrix& A ){
68 return gsl_linalg_householder_hm( tau, v.get(), A.get() ); }
76 inline int householder_mh( double tau, vector const& v, matrix& A ){
77 return gsl_linalg_householder_mh( tau, v.get(), A.get() ); }
85 inline int householder_hv( double tau, vector const& v, vector& w ){
86 return gsl_linalg_householder_hv( tau, v.get(), w.get() ); }
93 inline int householder_hm1( double tau, matrix& A ){
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() ); }
128 inline int hessenberg_decomp( matrix& A, vector& tau ){
129 return gsl_linalg_hessenberg_decomp( A.get(), tau.get() ); }
137 inline int hessenberg_unpack( matrix& H, vector& tau, matrix& U ){
138 return gsl_linalg_hessenberg_unpack( H.get(), tau.get(), U.get() ); }
146 inline int hessenberg_unpack_accum( matrix& H, vector& tau, matrix& U ){
147 return gsl_linalg_hessenberg_unpack_accum( H.get(), tau.get(), U.get() ); }
153 inline int hessenberg_set_zero( matrix& H ){ return gsl_linalg_hessenberg_set_zero( H.get() ); }
162 inline int hessenberg_submatrix( matrix& M, matrix& A, size_t top, vector& tau ){
163 return gsl_linalg_hessenberg_submatrix( M.get(), A.get(), top, tau.get() ); }
165
172 //inline int hessenberg( matrix& A, vector& tau ){ return gsl_linalg_hessenberg( A.get(), tau.get() ); }
182 inline int hesstri_decomp( matrix& A, matrix& B, matrix& U, matrix& V, vector& work ){
183 return gsl_linalg_hesstri_decomp( A.get(), B.get(), U.get(), V.get(), work.get() ); }
192 inline int SV_decomp( matrix& A, matrix& V, vector& S, vector& work ){
193 return gsl_linalg_SV_decomp( A.get(), V.get(), S.get(), work.get() ); }
203 inline int SV_decomp_mod( matrix& A, matrix& X, matrix& V, vector& S, vector& work ){
204 return gsl_linalg_SV_decomp_mod( A.get(), X.get(), V.get(), S.get(), work.get() ); }
212 inline int SV_decomp_jacobi( matrix& A, matrix& Q, vector& S ){
213 return gsl_linalg_SV_decomp_jacobi( A.get(), Q.get(), S.get() ); }
223 inline int SV_solve( matrix const& U, matrix const& Q, vector const& S, vector const& b, vector& x ){
224 return gsl_linalg_SV_solve( U.get(), Q.get(), S.get(), b.get(), x.get() ); }
225#ifndef DOXYGEN_SKIP
233 inline int LU_decomp( matrix& A, permutation& p, int* signum ){
234 return gsl_linalg_LU_decomp( A.get(), p.get(), signum ); }
235#endif
243 inline int LU_decomp( matrix& A, permutation& p, int& signum ){
244 return gsl_linalg_LU_decomp( A.get(), p.get(), &signum ); }
253 inline int LU_solve( matrix const& LU, permutation const& p, vector const& b, vector& x ){
254 return gsl_linalg_LU_solve( LU.get(), p.get(), b.get(), x.get() ); }
262 inline int LU_svx( matrix const& LU, permutation const& p, vector& x ){
263 return gsl_linalg_LU_svx( LU.get(), p.get(), x.get() ); }
274 inline int LU_refine( matrix const& A, matrix const& LU, permutation const& p, vector const& b,
275 vector& x, vector& residual ){
276 return gsl_linalg_LU_refine( A.get(), LU.get(), p.get(), b.get(), x.get(), residual.get() ); }
284 inline int LU_invert( matrix const& LU, permutation const& p, matrix& inverse ){
285 return gsl_linalg_LU_invert( LU.get(), p.get(), inverse.get() ); }
292 inline double LU_det( matrix& LU, int signum ){ return gsl_linalg_LU_det( LU.get(), signum ); }
299 inline double LU_lndet( matrix& LU ){ return gsl_linalg_LU_lndet( LU.get() ); }
306 inline int LU_sgndet( matrix& lu, int signum ){ return gsl_linalg_LU_sgndet( lu.get(), signum ); }
307#ifndef DOXYGEN_SKIP
315 inline int complex_LU_decomp( matrix_complex& A, permutation& p, int* signum ){
316 return gsl_linalg_complex_LU_decomp( A.get(), p.get(), signum ); }
317#endif
325 inline int complex_LU_decomp( matrix_complex& A, permutation& p, int& signum ){
326 return gsl_linalg_complex_LU_decomp( A.get(), p.get(), &signum ); }
335 inline int complex_LU_solve( matrix_complex const& LU, permutation const& p,
336 vector_complex const& b, vector_complex& x ){
337 return gsl_linalg_complex_LU_solve( LU.get(), p.get(), b.get(), x.get() ); }
345 inline int complex_LU_svx( matrix_complex const& LU, permutation const& p, vector_complex& x ){
346 return gsl_linalg_complex_LU_svx( LU.get(), p.get(), x.get() ); }
357 inline int complex_LU_refine( matrix_complex const& A, matrix_complex const& LU,
358 permutation const& p, vector_complex const& b,
360 return gsl_linalg_complex_LU_refine( A.get(), LU.get(), p.get(), b.get(), x.get(), residual.get() ); }
368 inline int complex_LU_invert( matrix_complex const& LU, permutation const& p,
370 return gsl_linalg_complex_LU_invert( LU.get(), p.get(), inverse.get() ); }
377 inline complex complex_LU_det( matrix_complex& LU, int signum ){
378 return gsl_linalg_complex_LU_det( LU.get(), signum ); }
384 inline double complex_LU_lndet( matrix_complex& LU ){
385 return gsl_linalg_complex_LU_lndet( LU.get() ); }
392 inline complex complex_LU_sgndet( matrix_complex& LU, int signum ){
393 return gsl_linalg_complex_LU_sgndet( LU.get(), signum ); }
400 inline int QR_decomp( matrix& A, vector& tau ){ return gsl_linalg_QR_decomp( A.get(), tau.get() ); }
409 inline int QR_solve( matrix const& QR, vector const& tau, vector const& b, vector& x ){
410 return gsl_linalg_QR_solve( QR.get(), tau.get(), b.get(), x.get() ); }
418 inline int QR_svx( matrix const& QR, vector const& tau, vector& x ){
419 return gsl_linalg_QR_svx( QR.get(), tau.get(), x.get() ); }
429 inline int QR_lssolve( matrix const& QR, vector const& tau, vector const& b, vector& x,
430 vector& residual ){
431 return gsl_linalg_QR_lssolve( QR.get(), tau.get(), b.get(), x.get(), residual.get() ); }
440 inline int QR_QRsolve( matrix& Q, matrix& R, vector const& b, vector& x ){
441 return gsl_linalg_QR_QRsolve( Q.get(), R.get(), b.get(), x.get() ); }
449 inline int QR_Rsolve( matrix const& QR, vector const& b, vector& x ){
450 return gsl_linalg_QR_Rsolve( QR.get(), b.get(), x.get() ); }
457 inline int QR_Rsvx( matrix const& QR, vector& x ){ return gsl_linalg_QR_Rsvx( QR.get(), x.get() ); }
466 inline int QR_update( matrix& Q, matrix& R, vector& w, vector const& v ){
467 return gsl_linalg_QR_update( Q.get(), R.get(), w.get(), v.get() ); }
475 inline int QR_QTvec( matrix const& QR, vector const& tau, vector& v ){
476 return gsl_linalg_QR_QTvec( QR.get(), tau.get(), v.get() ); }
484 inline int QR_Qvec( matrix const& QR, vector const& tau, vector& v ){
485 return gsl_linalg_QR_Qvec( QR.get(), tau.get(), v.get() ); }
493 inline int QR_QTmat( matrix const& QR, vector const& tau, matrix& A ){
494 return gsl_linalg_QR_QTmat( QR.get(), tau.get(), A.get() ); }
503 inline int QR_unpack( matrix const& QR, vector const& tau, matrix& Q, matrix& R ){
504 return gsl_linalg_QR_unpack( QR.get(), tau.get(), Q.get(), R.get() ); }
512 inline int R_solve( matrix const& R, vector const& b, vector& x ){
513 return gsl_linalg_R_solve( R.get(), b.get(), x.get() ); }
520 inline int R_svx( matrix const& R, vector& x ){
521 return gsl_linalg_R_svx( R.get(), x.get() ); }
522#ifndef DOXYGEN_SKIP
532 inline int QRPT_decomp( matrix& A, vector& tau, permutation& p, int* signum, vector& norm ){
533 return gsl_linalg_QRPT_decomp( A.get(), tau.get(), p.get(), signum, norm.get() ); }
534#endif
544 inline int QRPT_decomp( matrix& A, vector& tau, permutation& p, int& signum, vector& norm ){
545 return gsl_linalg_QRPT_decomp( A.get(), tau.get(), p.get(), &signum, norm.get() ); }
546#ifndef DOXYGEN_SKIP
558 inline int QRPT_decomp2( matrix const& A, matrix& q, matrix& r, vector& tau,
559 permutation& p, int* signum, vector& norm ){
560 return gsl_linalg_QRPT_decomp2( A.get(), q.get(), r.get(), tau.get(), p.get(), signum, norm.get() ); }
561#endif
573 inline int QRPT_decomp2( matrix const& A, matrix& q, matrix& r, vector& tau,
574 permutation& p, int& signum, vector& norm ){
575 return gsl_linalg_QRPT_decomp2( A.get(), q.get(), r.get(), tau.get(), p.get(), &signum, norm.get() ); }
585 inline int QRPT_solve( matrix const& QR, vector const& tau, permutation const& p,
586 vector const& b, vector& x ){
587 return gsl_linalg_QRPT_solve( QR.get(), tau.get(), p.get(), b.get(), x.get() ); }
596 inline int QRPT_svx( matrix const& QR, vector const& tau, permutation const& p, vector& x ){
597 return gsl_linalg_QRPT_svx( QR.get(), tau.get(), p.get(), x.get() ); }
607 inline int QRPT_QRsolve( matrix const& Q, matrix const& R, permutation const& p,
608 vector const& b, vector& x ){
609 return gsl_linalg_QRPT_QRsolve( Q.get(), R.get(), p.get(), b.get(), x.get() ); }
618 inline int QRPT_Rsolve( matrix const& QR, permutation const& p, vector const& b, vector& x ){
619 return gsl_linalg_QRPT_Rsolve( QR.get(), p.get(), b.get(), x.get() ); }
627 inline int QRPT_Rsvx( matrix const& QR, permutation const& p, vector& x ){
628 return gsl_linalg_QRPT_Rsvx( QR.get(), p.get(), x.get() ); }
638 inline int QRPT_update( matrix& Q, matrix& R, permutation const& p, vector& u, vector const& v ){
639 return gsl_linalg_QRPT_update( Q.get(), R.get(), p.get(), u.get(), v.get() ); }
646 inline int LQ_decomp( matrix& A, vector& tau ){
647 return gsl_linalg_LQ_decomp( A.get(), tau.get() ); }
656 inline int LQ_solve_T( matrix const& LQ, vector const& tau, vector const& b, vector& x ){
657 return gsl_linalg_LQ_solve_T( LQ.get(), tau.get(), b.get(), x.get() ); }
665 inline int LQ_svx_T( matrix const& LQ, vector const& tau, vector& x ){
666 return gsl_linalg_LQ_svx_T( LQ.get(), tau.get(), x.get() ); }
676 inline int LQ_lssolve_T( matrix const& LQ, vector const& tau, vector const& b,
677 vector& x, vector& residual ){
678 return gsl_linalg_LQ_lssolve_T( LQ.get(), tau.get(), b.get(), x.get(), residual.get() ); }
686 inline int LQ_Lsolve_T( matrix const& LQ, vector const& b, vector& x ){
687 return gsl_linalg_LQ_Lsolve_T( LQ.get(), b.get(), x.get() ); }
694 inline int LQ_Lsvx_T( matrix const& LQ, vector& x ){
695 return gsl_linalg_LQ_Lsvx_T( LQ.get(), x.get() ); }
703 inline int L_solve_T( matrix const& L, vector const& b, vector& x ){
704 return gsl_linalg_L_solve_T( L.get(), b.get(), x.get() ); }
712 inline int LQ_vecQ( matrix const& LQ, vector const& tau, vector& v ){
713 return gsl_linalg_LQ_vecQ( LQ.get(), tau.get(), v.get() ); }
721 inline int LQ_vecQT( matrix const& LQ, vector const& tau, vector& v ){
722 return gsl_linalg_LQ_vecQT( LQ.get(), tau.get(), v.get() ); }
731 inline int LQ_unpack( matrix const& LQ, vector const& tau, matrix& Q, matrix& L ){
732 return gsl_linalg_LQ_unpack( LQ.get(), tau.get(), Q.get(), L.get() ); }
741 inline int LQ_update( matrix& Q, matrix& R, vector const& v, vector& w ){
742 return gsl_linalg_LQ_update( Q.get(), R.get(), v.get(), w.get() ); }
751 inline int LQ_LQsolve( matrix& Q, matrix& L, vector const& b, vector& x ){
752 return gsl_linalg_LQ_LQsolve( Q.get(), L.get(), b.get(), x.get() ); }
753#ifndef DOXYGEN_SKIP
763 inline int PTLQ_decomp( matrix& A, vector& tau, permutation& p, int* signum, vector& norm ){
764 return gsl_linalg_PTLQ_decomp( A.get(), tau.get(), p.get(), signum, norm.get() ); }
765#endif
775 inline int PTLQ_decomp( matrix& A, vector& tau, permutation& p, int& signum, vector& norm ){
776 return gsl_linalg_PTLQ_decomp( A.get(), tau.get(), p.get(), &signum, norm.get() ); }
788 inline int PTLQ_decomp2( matrix const& A, matrix& q, matrix& r, vector& tau,
789 permutation& p, int* signum, vector& norm ){
790 return gsl_linalg_PTLQ_decomp2( A.get(), q.get(), r.get(), tau.get(), p.get(), signum, norm.get() ); }
802 inline int PTLQ_decomp2( matrix const& A, matrix& q, matrix& r, vector& tau,
803 permutation& p, int& signum, vector& norm ){
804 return gsl_linalg_PTLQ_decomp2( A.get(), q.get(), r.get(), tau.get(), p.get(), &signum, norm.get() ); }
814 inline int PTLQ_solve_T( matrix const& QR, vector const& tau, permutation const& p,
815 vector const& b, vector& x ){
816 return gsl_linalg_PTLQ_solve_T( QR.get(), tau.get(), p.get(), b.get(), x.get() ); }
825 inline int PTLQ_svx_T( matrix const& LQ, vector const& tau, permutation const& p, vector& x ){
826 return gsl_linalg_PTLQ_svx_T( LQ.get(), tau.get(), p.get(), x.get() ); }
836 inline int PTLQ_LQsolve_T( matrix const& Q, matrix const& L, permutation const& p,
837 vector const& b, vector& x ){
838 return gsl_linalg_PTLQ_LQsolve_T( Q.get(), L.get(), p.get(), b.get(), x.get() ); }
847 inline int PTLQ_Lsolve_T( matrix const& LQ, permutation const& p, vector const& b, vector& x ){
848 return gsl_linalg_PTLQ_Lsolve_T( LQ.get(), p.get(), b.get(), x.get() ); }
856 inline int PTLQ_Lsvx_T( matrix const& LQ, permutation const& p, vector& x ){
857 return gsl_linalg_PTLQ_Lsvx_T( LQ.get(), p.get(), x.get() ); }
867 inline int PTLQ_update( matrix& Q, matrix& L, permutation const& p, vector const& v, vector& w ){
868 return gsl_linalg_PTLQ_update( Q.get(), L.get(), p.get(), v.get(), w.get() ); }
874 inline int cholesky_decomp( matrix& A ){ return gsl_linalg_cholesky_decomp( A.get() ); }
882 inline int cholesky_solve( matrix const& cholesky, vector const& b, vector& x ){
883 return gsl_linalg_cholesky_solve( cholesky.get(), b.get(), x.get() ); }
890 inline int cholesky_svx( matrix const& cholesky, vector& x ){
891 return gsl_linalg_cholesky_svx( cholesky.get(), x.get() ); }
897 inline int cholesky_invert( matrix& cholesky ){ return gsl_linalg_cholesky_invert( cholesky.get() ); }
905 return gsl_linalg_cholesky_decomp_unit( A.get(), D.get() ); }
912 return gsl_linalg_complex_cholesky_decomp( A.get() ); }
920 inline int complex_cholesky_solve( matrix_complex const& cholesky, vector_complex const& b,
921 vector_complex& x ){
922 return gsl_linalg_complex_cholesky_solve( cholesky.get(), b.get(), x.get() ); }
929 inline int complex_cholesky_svx( matrix_complex const& cholesky, vector_complex& x ){
930 return gsl_linalg_complex_cholesky_svx( cholesky.get(), x.get() ); }
937 inline int symmtd_decomp( matrix& A, vector& tau ){
938 return gsl_linalg_symmtd_decomp( A.get(), tau.get() ); }
948 inline int symmtd_unpack( matrix const& A, vector const& tau, matrix& Q,
949 vector& diag, vector& subdiag ){
950 return gsl_linalg_symmtd_unpack( A.get(), tau.get(), Q.get(), diag.get(), subdiag.get() ); }
958 inline int symmtd_unpack_T( matrix const& A, vector& diag, vector& subdiag ){
959 return gsl_linalg_symmtd_unpack_T( A.get(), diag.get(), subdiag.get() ); }
967 return gsl_linalg_hermtd_decomp( A.get(), tau.get() ); }
977 inline int hermtd_unpack( matrix_complex const& A, vector_complex const& tau,
978 matrix_complex& U, vector& diag, vector& sudiag ){
979 return gsl_linalg_hermtd_unpack( A.get(), tau.get(), U.get(), diag.get(), sudiag.get() ); }
987 inline int hermtd_unpack_T( matrix_complex const& A, vector& diag, vector& subdiag ){
988 return gsl_linalg_hermtd_unpack_T( A.get(), diag.get(), subdiag.get() ); }
996 inline int HH_solve( matrix& A, vector const& b, vector& x ){
997 return gsl_linalg_HH_solve( A.get(), b.get(), x.get() ); }
1004 inline int HH_svx( matrix& A, vector& x ){ return gsl_linalg_HH_svx( A.get(), x.get() ); }
1013 inline int solve_symm_tridiag( vector const& diag, vector const& offdiag, vector const& b,
1014 vector& x ){
1015 return gsl_linalg_solve_symm_tridiag( diag.get(), offdiag.get(), b.get(), x.get() ); }
1025 inline int solve_tridiag( vector const& diag, vector const& abovediag, vector const& belowdiag,
1026 vector const& b, vector& x ){
1027 return gsl_linalg_solve_tridiag( diag.get(), abovediag.get(), belowdiag.get(), b.get(), x.get() ); }
1036 inline int solve_symm_cyc_tridiag( vector const& diag, vector const& offdiag,
1037 vector const& b, vector& x ){
1038 return gsl_linalg_solve_symm_cyc_tridiag( diag.get(), offdiag.get(), b.get(), x.get() ); }
1048 inline int solve_cyc_tridiag( vector const& diag, vector const& abovediag, vector const& belowdiag,
1049 vector const& b, vector& x ){
1050 return gsl_linalg_solve_cyc_tridiag( diag.get(), abovediag.get(), belowdiag.get(), b.get(), x.get() ); }
1058 inline int bidiag_decomp( matrix& A, vector& tau_U, vector& tau_V ){
1059 return gsl_linalg_bidiag_decomp( A.get(), tau_U.get(), tau_V.get() ); }
1071 inline int bidiag_unpack( matrix const& A, vector const& tau_U, matrix& U,
1072 vector const& tau_V, matrix& V, vector& diag, vector& superdiag ){
1073 return gsl_linalg_bidiag_unpack( A.get(), tau_U.get(), U.get(), tau_V.get(),
1074 V.get(), diag.get(), superdiag.get() ); }
1083 inline int bidiag_unpack2( matrix& A, vector& tau_U, vector& tau_V, matrix& V ){
1084 return gsl_linalg_bidiag_unpack2( A.get(), tau_U.get(), tau_V.get(), V.get() ); }
1092 inline int bidiag_unpack_B( matrix const& A, vector& diag, vector& superdiag ){
1093 return gsl_linalg_bidiag_unpack_B( A.get(), diag.get(), superdiag.get() ); }
1100 inline int balance_matrix( matrix& A, vector& D ){
1101 return gsl_linalg_balance_matrix( A.get(), D.get() ); }
1108 inline int balance_accum( matrix& A, vector& D ){
1109 return gsl_linalg_balance_accum( A.get(), D.get() ); }
1116 inline int balance_columns( matrix& A, vector& D ){
1117 return gsl_linalg_balance_columns( A.get(), D.get() ); }
1124 inline int SV_leverage( matrix const& U, vector& h){
1125 return gsl_linalg_SV_leverage( U.get(), h.get() ); }
1134 inline int householder_left( double const tau, vector const& v,
1135 matrix& A, vector& work ){
1136 return gsl_linalg_householder_left( tau, v.get(), A.get(), work.get() ); }
1145 inline int householder_right( double const tau, vector const& v, matrix& A, vector& work ){
1146 return gsl_linalg_householder_right( tau, v.get(), A.get(), work.get() ); }
1155 inline int complex_householder_left( complex const tau, vector_complex const& v,
1156 matrix_complex& A, vector_complex& work ){
1157 return gsl_linalg_complex_householder_left( tau.get(), v.get(), A.get(), work.get() ); }
1164 inline int LU_invx( matrix& LU, permutation const& p ){
1165 return gsl_linalg_LU_invx( LU.get(), p.get() ); }
1172 inline int complex_LU_invx( matrix_complex& LU, permutation const& p ){
1173 return gsl_linalg_complex_LU_invx( LU.get(), p.get() ); }
1182 inline int QR_QTmat_r( matrix const& QR, matrix const& T, matrix& B, matrix& work ){
1183 return gsl_linalg_QR_QTmat_r( QR.get(), T.get(), B.get(), work.get() ); }
1191 inline int QR_matQ( matrix const& QR, vector const& tau, matrix& A ){
1192 return gsl_linalg_QR_matQ( QR.get(), tau.get(), A.get() ); }
1201 inline int QR_unpack_r( matrix const& QR, matrix const& T, matrix& Q, matrix& R ){
1202 return gsl_linalg_QR_unpack_r( QR.get(), T.get(), Q.get(), R.get() ); }
1210 inline int QR_rcond( matrix const& QR, double& rcond, vector& work ){
1211 return gsl_linalg_QR_rcond( QR.get(), &rcond, work.get() ); }
1222 inline int QRPT_lssolve( matrix const& QR, vector const& tau, permutation const& p,
1223 vector const& b, vector& x, vector& residual ){
1224 return gsl_linalg_QRPT_lssolve( QR.get(), tau.get(), p.get(), b.get(), x.get(),
1225 residual.get() ); }
1237 inline int QRPT_lssolve2( matrix const& QR, vector const& tau, permutation const& p,
1238 vector const& b, size_t const rank, vector& x, vector& residual ){
1239 return gsl_linalg_QRPT_lssolve2( QR.get(), tau.get(), p.get(), b.get(), rank, x.get(),
1240 residual.get() ); }
1247 inline size_t QRPT_rank( matrix const& QR, double const tol ){
1248 return gsl_linalg_QRPT_rank( QR.get(), tol ); }
1256 inline int QRPT_rcond( matrix const& QR, double& rcond, vector& work ){
1257 return gsl_linalg_QRPT_rcond( QR.get(), &rcond, work.get() ); }
1268 inline int COD_decomp( matrix& A, vector& tau_Q, vector& tau_Z, permutation& p,
1269 size_t& rank, vector& work ){
1270 return gsl_linalg_COD_decomp( A.get(), tau_Q.get(), tau_Z.get(), p.get(), &rank,
1271 work.get() ); }
1283 inline int COD_decomp_e( matrix& A, vector& tau_Q, vector& tau_Z, permutation& p,
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() ); }
1299 inline int COD_lssolve( matrix const& QRZT, vector const& tau_Q, vector const& tau_Z,
1300 permutation const& perm, size_t const rank, vector const& b,
1301 vector& x, vector& residual ){
1302 return gsl_linalg_COD_lssolve( QRZT.get(), tau_Q.get(), tau_Z.get(), perm.get(),
1303 rank, b.get(), x.get(), residual.get() ); }
1319 inline int COD_lssolve2( double const lambda, matrix const& QRZT, vector const& tau_Q,
1320 vector const& tau_Z, permutation const& perm, size_t const rank,
1321 vector const& b, vector& x, vector& residual, matrix& S,
1322 vector& work ){
1323 return gsl_linalg_COD_lssolve2( lambda, QRZT.get(), tau_Q.get(), tau_Z.get(), perm.get(),
1324 rank, b.get(), x.get(), residual.get(), S.get(),
1325 work.get() ); }
1337 inline int COD_unpack( matrix const& QRZT, vector const& tau_Q, vector const& tau_Z,
1338 size_t const rank, matrix& Q, matrix& R, matrix& Z ){
1339 return gsl_linalg_COD_unpack( QRZT.get(), tau_Q.get(), tau_Z.get(), rank, Q.get(),
1340 R.get(), Z.get() ); }
1350 inline int COD_matZ( matrix const& QRZT, vector const& tau_Z, size_t const rank,
1351 matrix& A, vector& work ){
1352 return gsl_linalg_COD_matZ( QRZT.get(), tau_Z.get(), rank, A.get(), work.get() ); }
1362 inline int LQ_lssolve( matrix const& LQ, vector const& tau, vector const& b, vector& x,
1363 vector& residual ){
1364 return gsl_linalg_LQ_lssolve( LQ.get(), tau.get(), b.get(), x.get(), residual.get() ); }
1372 inline int LQ_QTvec( matrix const& LQ, vector const& tau, vector& v ){
1373 return gsl_linalg_LQ_QTvec( LQ.get(), tau.get(), v.get() ); }
1379 inline int cholesky_decomp1( matrix& A ){ return gsl_linalg_cholesky_decomp1( A.get() ); }
1387 inline int cholesky_solve_mat( matrix const& cholesky, matrix const& B, matrix& X ){
1388 return gsl_linalg_cholesky_solve_mat( cholesky.get(), B.get(), X.get() ); }
1395 inline int cholesky_svx_mat( matrix const& cholesky, matrix& X ){
1396 return gsl_linalg_cholesky_svx_mat( cholesky.get(), X.get() ); }
1403 inline int cholesky_scale( matrix const& A, vector& S ){
1404 return gsl_linalg_cholesky_scale( A.get(), S.get() ); }
1411 inline int cholesky_scale_apply( matrix& A, vector const& S ){
1412 return gsl_linalg_cholesky_scale_apply( A.get(), S.get() ); }
1419 inline int cholesky_decomp2( matrix& A, vector& S ){
1420 return gsl_linalg_cholesky_decomp2( A.get(), S.get() ); }
1428 inline int cholesky_svx2( matrix const& LLT, vector const& S, vector& x ){
1429 return gsl_linalg_cholesky_svx2( LLT.get(), S.get(), x.get() ); }
1438 inline int cholesky_solve2( matrix const& LLT, vector const& S, vector const& b,
1439 vector& x ){
1440 return gsl_linalg_cholesky_solve2( LLT.get(), S.get(), b.get(), x.get() ); }
1448 inline int cholesky_rcond( matrix const& LLT, double& rcond, vector& work ){
1449 return gsl_linalg_cholesky_rcond( LLT.get(), &rcond, work.get() ); }
1457 return gsl_linalg_pcholesky_decomp( A.get(), p.get() ); }
1466 inline int pcholesky_solve( matrix const& LDLT, permutation const& p, vector const& b,
1467 vector& x ){
1468 return gsl_linalg_pcholesky_solve( LDLT.get(), p.get(), b.get(), x.get() ); }
1476 inline int pcholesky_svx( matrix const& LDLT, permutation const& p, vector& x ){
1477 return gsl_linalg_pcholesky_svx( LDLT.get(), p.get(), x.get() ); }
1485 inline int pcholesky_decomp2( matrix& A, permutation& p, vector& S ){
1486 return gsl_linalg_pcholesky_decomp2( A.get(), p.get(), S.get() ); }
1496 inline int pcholesky_solve2( matrix const& LDLT, permutation const& p, vector const& S,
1497 vector const& b, vector& x ){
1498 return gsl_linalg_pcholesky_solve2( LDLT.get(), p.get(), S.get(), b.get(), x.get() ); }
1507 inline int pcholesky_svx2( matrix const& LDLT, permutation const& p, vector const& S,
1508 vector& x ){
1509 return gsl_linalg_pcholesky_svx2( LDLT.get(), p.get(), S.get(), x.get() ); }
1517 inline int pcholesky_invert( matrix const& LDLT, permutation const& p, matrix& Ainv ){
1518 return gsl_linalg_pcholesky_invert( LDLT.get(), p.get(), Ainv.get() ); }
1527 inline int pcholesky_rcond( matrix const& LDLT, permutation const& p, double& rcond,
1528 vector& work ){
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() ); }
1547 inline int mcholesky_solve( matrix const& LDLT, permutation const& p, vector const& b,
1548 vector& x ){
1549 return gsl_linalg_mcholesky_solve( LDLT.get(), p.get(), b.get(), x.get() ); }
1557 inline int mcholesky_svx( matrix const& LDLT, permutation const& p, vector& x ){
1558 return gsl_linalg_mcholesky_svx( LDLT.get(), p.get(), x.get() ); }
1567 inline int mcholesky_rcond( matrix const& LDLT, permutation const& p, double& rcond,
1568 vector& work ){
1569 return gsl_linalg_mcholesky_rcond( LDLT.get(), p.get(), &rcond, work.get() ); }
1577 inline int mcholesky_invert( matrix const& LDLT, permutation const& p, matrix& Ainv ){
1578 return gsl_linalg_mcholesky_invert( LDLT.get(), p.get(), Ainv.get() ); }
1585 inline int
1587 return gsl_linalg_cholesky_band_decomp( A.get() ); }
1595 inline int cholesky_band_solve( matrix const& LLT, vector const& b, vector& x ){
1596 return gsl_linalg_cholesky_band_solve( LLT.get(), b.get(), x.get() ); }
1603 inline int cholesky_band_svx( matrix const& LLT, vector& x ){
1604 return gsl_linalg_cholesky_band_svx( LLT.get(), x.get() ); }
1611 inline int cholesky_band_invert( matrix const& LLT, matrix& Ainv ){
1612 return gsl_linalg_cholesky_band_invert( LLT.get(), Ainv.get() ); }
1619 inline int cholesky_band_unpack( matrix const& LLT, matrix& L ){
1620 return gsl_linalg_cholesky_band_unpack( LLT.get(), L.get() ); }
1628 inline int cholesky_band_rcond( matrix const& LLT, double& rcond, vector& work ){
1629 return gsl_linalg_cholesky_band_rcond( LLT.get(), &rcond, work.get() ); }
1635 inline int ldlt_decomp( matrix& A ){ return gsl_linalg_ldlt_decomp( A.get() ); }
1643 inline int ldlt_solve( matrix const& LDLT, vector const& b, vector& x ){
1644 return gsl_linalg_ldlt_solve( LDLT.get(), b.get(), x.get() ); }
1651 inline int ldlt_svx( matrix const& LDLT, vector& x ){
1652 return gsl_linalg_ldlt_svx( LDLT.get(), x.get() ); }
1660 inline int ldlt_rcond( matrix const& LDLT, double& rcond, vector& work ){
1661 return gsl_linalg_ldlt_rcond( LDLT.get(), &rcond, work.get() ); }
1668 inline int ldlt_band_decomp( matrix& A ){
1669 return gsl_linalg_ldlt_band_decomp( A.get() ); }
1677 inline int ldlt_band_solve( matrix const& LDLT, vector const& b, vector& x ){
1678 return gsl_linalg_ldlt_band_solve( LDLT.get(), b.get(), x.get() ); }
1685 inline int ldlt_band_svx( matrix const& LDLT, vector& x ){
1686 return gsl_linalg_ldlt_band_svx( LDLT.get(), x.get() ); }
1694 inline int ldlt_band_unpack( matrix const& LDLT, matrix& L, vector& D ){
1695 return gsl_linalg_ldlt_band_unpack( LDLT.get(), L.get(), D.get() ); }
1703 inline int ldlt_band_rcond( matrix const& LDLT, double* rcond, vector& work ){
1704 return gsl_linalg_ldlt_band_rcond( LDLT.get(), rcond, work.get() ); }
1712 inline int ldlt_band_rcond( matrix const& LDLT, double& rcond, vector& work ){
1713 return gsl_linalg_ldlt_band_rcond( LDLT.get(), &rcond, work.get() ); }
1722 inline int tri_rcond( CBLAS_UPLO_t Uplo, matrix const& A, double& rcond, vector& work ){
1723 return gsl_linalg_tri_rcond( Uplo, A.get(), &rcond, work.get() ); }
1731 inline int tri_upper_rcond( matrix const& A, double& rcond, vector& work ){
1732 return gsl_linalg_tri_upper_rcond( A.get(), &rcond, work.get() ); }
1740 inline int tri_lower_rcond( matrix const& A, double& rcond, vector& work ){
1741 return gsl_linalg_tri_lower_rcond( A.get(), &rcond, work.get() ); }
1742#ifndef DOXYGEN_SKIP
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() ); }
1756#endif
1762 inline int tri_upper_invert( matrix& T ){ return gsl_linalg_tri_upper_invert( T.get() ); }
1768 inline int
1769 tri_lower_invert( matrix& T ){ return gsl_linalg_tri_lower_invert( T.get() ); }
1776 return gsl_linalg_tri_upper_unit_invert( T.get() ); }
1783 return gsl_linalg_tri_lower_unit_invert( T.get() ); }
1791 inline int tri_invert( CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix& T ){
1792 return gsl_linalg_tri_invert( Uplo, Diag, T.get() ); }
1800 inline int complex_tri_invert( CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix_complex& T ){
1801 return gsl_linalg_complex_tri_invert( Uplo, Diag, T.get() ); }
1807 inline int tri_LTL( matrix& L ){ return gsl_linalg_tri_LTL( L.get() ); }
1813 inline int tri_UL( matrix& LU ){ return gsl_linalg_tri_UL( LU.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 ); }
1845 inline void givens_gv( vector& v, size_t const i, size_t const j, double const c,
1846 double const s ){
1847 return gsl_linalg_givens_gv( v.get(), i, j, c, s ); }
1857 inline int LU_band_decomp( size_t const M, size_t const lb, size_t const ub,
1858 matrix& AB, vector_uint& piv ){
1859 return gsl_linalg_LU_band_decomp( M, lb, ub, AB.get(), piv.get() ); }
1870 inline int LU_band_solve( size_t const lb, size_t const ub,
1871 matrix const& LUB, vector_uint const& piv,
1872 vector const& b, vector& x ){
1873 return gsl_linalg_LU_band_solve( lb, ub, LUB.get(), piv.get(), b.get(), x.get() ); }
1883 inline int LU_band_svx( size_t const lb, size_t const ub, matrix const& LUB,
1884 vector_uint const& piv, vector& x ){
1885 return gsl_linalg_LU_band_svx( lb, ub, LUB.get(), piv.get(), x.get() ); }
1897 inline int LU_band_unpack( size_t const M, size_t const lb, size_t const ub,
1898 matrix const& LUB, vector_uint const& piv,
1899 matrix& L, matrix& U ){
1900 return gsl_linalg_LU_band_unpack( M, lb, ub, LUB.get(), piv.get(), L.get(), U.get() ); }
1907 inline int QR_decomp_old( matrix& A, vector& tau ){
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() ); }
1933 inline int complex_QR_solve( matrix_complex const& QR, vector_complex const& tau,
1934 vector_complex const& b, vector_complex& x ){
1935 return gsl_linalg_complex_QR_solve( QR.get(), tau.get(), b.get(), x.get() ); }
1944 inline int complex_QR_solve_r( matrix_complex const& QR, matrix_complex const& T,
1945 vector_complex const& b, vector_complex& x ){
1946 return gsl_linalg_complex_QR_solve_r( QR.get(), T.get(), b.get(), x.get() ); }
1954 inline int complex_QR_svx( matrix_complex const& QR, vector_complex const& tau,
1955 vector_complex& x ){
1956 return gsl_linalg_complex_QR_svx( QR.get(), tau.get(), x.get() ); }
1966 inline int complex_QR_lssolve( matrix_complex const& QR, vector_complex const& tau,
1967 vector_complex const& b, vector_complex& x,
1969 return gsl_linalg_complex_QR_lssolve( QR.get(), tau.get(), b.get(), x.get(),
1970 residual.get() ); }
1980 inline int complex_QR_lssolve_r( matrix_complex const& QR, matrix_complex const& T,
1981 vector_complex const& b, vector_complex& x,
1982 vector_complex& work ){
1983 return gsl_linalg_complex_QR_lssolve_r( QR.get(), T.get(), b.get(), x.get(),
1984 work.get() ); }
1992 inline int complex_QR_QHvec( matrix_complex const& QR, vector_complex const& tau,
1993 vector_complex& v ){
1994 return gsl_linalg_complex_QR_QHvec( QR.get(), tau.get(), v.get() ); }
2003 inline int complex_QR_QHvec_r( matrix_complex const& QR, matrix_complex const& T,
2005 return gsl_linalg_complex_QR_QHvec_r( QR.get(), T.get(), b.get(), work.get() ); }
2013 inline int complex_QR_Qvec( matrix_complex const& QR, vector_complex const& tau,
2014 vector_complex& v ){
2015 return gsl_linalg_complex_QR_Qvec( QR.get(), tau.get(), v.get() ); }
2024 inline int complex_QR_unpack( matrix_complex const& QR, vector_complex const& tau,
2026 return gsl_linalg_complex_QR_unpack( QR.get(), tau.get(), Q.get(), R.get() ); }
2035 inline int complex_QR_unpack_r( matrix_complex const& QR, matrix_complex const& T,
2037 return gsl_linalg_complex_QR_unpack_r( QR.get(), T.get(), Q.get(), R.get() ); }
2047 inline int QR_band_decomp_L2( size_t const M, size_t const p, size_t const q, matrix& AB,
2048 vector& tau ){
2049 return gsl_linalg_QR_band_decomp_L2( M, p, q, AB.get(), tau.get() ); }
2060 inline int QR_band_unpack_L2( size_t const p, size_t const q, matrix const& QRB,
2061 vector const& tau, matrix& Q, matrix& R ){
2062 return gsl_linalg_QR_band_unpack_L2( p, q, QRB.get(), tau.get(), Q.get(), R.get() ); }
2071 inline int QR_UD_decomp( matrix& U, vector const& D, matrix& Y, matrix& T ){
2072 return gsl_linalg_QR_UD_decomp( U.get(), D.get(), Y.get(), T.get() ); }
2083 inline int QR_UD_lssolve( matrix const& R, matrix const& Y, matrix const& T, vector const& b,
2084 vector& x, vector& work ){
2085 return gsl_linalg_QR_UD_lssolve( R.get(), Y.get(), T.get(), b.get(), x.get(),
2086 work.get() ); }
2094 inline int QR_UR_decomp( matrix& S, matrix& A, matrix& T ){
2095 return gsl_linalg_QR_UR_decomp( S.get(), A.get(), T.get() ); }
2103 inline int QR_UU_decomp( matrix& U, matrix& S, matrix& T ){
2104 return gsl_linalg_QR_UU_decomp( U.get(), S.get(), T.get() ); }
2115 inline int QR_UU_lssolve( matrix const& R, matrix const& Y, matrix const& T, vector const& b,
2116 vector& x, vector& work ){
2117 return gsl_linalg_QR_UU_lssolve( R.get(), Y.get(), T.get(), b.get(), x.get(),
2118 work.get() ); }
2127 inline int QR_UU_QTvec( matrix const& Y, matrix const& T, vector& b, vector& work ){
2128 return gsl_linalg_QR_UU_QTvec( Y.get(), T.get(), b.get(), work.get() ); }
2136 inline int QR_UZ_decomp( matrix& S, matrix& A, matrix& T ){
2137 return gsl_linalg_QR_UZ_decomp( S.get(), A.get(), T.get() ); }
2144 inline int QL_decomp( matrix& A, vector& tau ){
2145 return gsl_linalg_QL_decomp( A.get(), tau.get() ); }
2154 inline int QL_unpack( matrix const& QL, vector const& tau, matrix& Q, matrix& L ){
2155 return gsl_linalg_QL_unpack( QL.get(), tau.get(), Q.get(), L.get() ); }
2156 }
2157}
2158
2159#endif
This class handles complex numbers.
Definition: complex.hpp:42
gsl_complex & get()
Get the base class object.
Definition: complex.hpp:106
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.
Definition: matrix.hpp:72
gsl_matrix * get()
Get the gsl_matrix.
Definition: matrix.hpp:1207
This class handles GSL permutation objects.
Definition: permutation.hpp:33
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.
Definition: vector_uint.hpp:45
gsl_vector_uint * get()
Get the gsl_vector_uint.
This class handles vector objects as shared handles.
Definition: vector.hpp:74
gsl_vector * get()
Get the gsl_vector.
Definition: vector.hpp:1320
complex inverse(complex const &a)
C++ version of gsl_complex_inverse().
Definition: complex.hpp:919
double LU_lndet(matrix &LU)
C++ version of gsl_linalg_LU_lndet().
Definition: linalg.hpp:299
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().
Definition: linalg.hpp:1222
int cholesky_svx2(matrix const &LLT, vector const &S, vector &x)
C++ version of gsl_linalg_cholesky_svx2().
Definition: linalg.hpp:1428
int mcholesky_rcond(matrix const &LDLT, permutation const &p, double &rcond, vector &work)
C++ version of gsl_linalg_mcholesky_rcond().
Definition: linalg.hpp:1567
int LU_solve(matrix const &LU, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_LU_solve().
Definition: linalg.hpp:253
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().
Definition: linalg.hpp:2024
int pcholesky_decomp(matrix &A, permutation &p)
C++ version of gsl_linalg_pcholesky_decomp().
Definition: linalg.hpp:1456
int cholesky_scale_apply(matrix &A, vector const &S)
C++ version of gsl_linalg_cholesky_scale_apply().
Definition: linalg.hpp:1411
int balance_matrix(matrix &A, vector &D)
C++ version of gsl_linalg_balance_matrix().
Definition: linalg.hpp:1100
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().
Definition: linalg.hpp:1883
int bidiag_decomp(matrix &A, vector &tau_U, vector &tau_V)
C++ version of gsl_linalg_bidiag_decomp().
Definition: linalg.hpp:1058
int QR_Qvec(matrix const &QR, vector const &tau, vector &v)
C++ version of gsl_linalg_QR_Qvec().
Definition: linalg.hpp:484
int cholesky_band_invert(matrix const &LLT, matrix &Ainv)
C++ version of gsl_linalg_cholesky_band_invert().
Definition: linalg.hpp:1611
int complex_QR_decomp(matrix_complex &A, vector_complex &tau)
C++ version of gsl_linalg_complex_QR_decomp().
Definition: linalg.hpp:1915
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().
Definition: linalg.hpp:2003
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().
Definition: linalg.hpp:2035
int cholesky_solve(matrix const &cholesky, vector const &b, vector &x)
C++ version of gsl_linalg_cholesky_solve().
Definition: linalg.hpp:882
int ldlt_svx(matrix const &LDLT, vector &x)
C++ version of gsl_linalg_ldlt_svx().
Definition: linalg.hpp:1651
int LU_decomp(matrix &A, permutation &p, int &signum)
C++ version of gsl_linalg_LU_decomp().
Definition: linalg.hpp:243
int balance_columns(matrix &A, vector &D)
C++ version of gsl_linalg_balance_columns().
Definition: linalg.hpp:1116
int R_svx(matrix const &R, vector &x)
C++ version of gsl_linalg_R_svx().
Definition: linalg.hpp:520
int pcholesky_solve(matrix const &LDLT, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_pcholesky_solve().
Definition: linalg.hpp:1466
int QRPT_svx(matrix const &QR, vector const &tau, permutation const &p, vector &x)
C++ version of gsl_linalg_QRPT_svx().
Definition: linalg.hpp:596
int solve_symm_tridiag(vector const &diag, vector const &offdiag, vector const &b, vector &x)
C++ version of gsl_linalg_solve_symm_tridiag().
Definition: linalg.hpp:1013
complex complex_LU_det(matrix_complex &LU, int signum)
C++ version of gsl_linalg_complex_LU_det().
Definition: linalg.hpp:377
int cholesky_svx_mat(matrix const &cholesky, matrix &X)
C++ version of gsl_linalg_cholesky_svx_mat().
Definition: linalg.hpp:1395
int cholesky_band_unpack(matrix const &LLT, matrix &L)
C++ version of gsl_linalg_cholesky_band_unpack().
Definition: linalg.hpp:1619
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().
Definition: linalg.hpp:1299
complex complex_LU_sgndet(matrix_complex &LU, int signum)
C++ version of gsl_linalg_complex_LU_sgndet().
Definition: linalg.hpp:392
int ldlt_band_rcond(matrix const &LDLT, double *rcond, vector &work)
C++ version of gsl_linalg_ldlt_band_rcond().
Definition: linalg.hpp:1703
int hessenberg_set_zero(matrix &H)
C++ version of gsl_linalg_hessenberg_set_zero().
Definition: linalg.hpp:153
int mcholesky_svx(matrix const &LDLT, permutation const &p, vector &x)
C++ version of gsl_linalg_mcholesky_svx().
Definition: linalg.hpp:1557
int tri_invert(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix &T)
C++ version of gsl_linalg_tri_invert().
Definition: linalg.hpp:1791
int complex_QR_QHvec(matrix_complex const &QR, vector_complex const &tau, vector_complex &v)
C++ version of gsl_linalg_complex_QR_QHvec().
Definition: linalg.hpp:1992
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().
Definition: linalg.hpp:1283
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().
Definition: linalg.hpp:1980
int QR_Rsolve(matrix const &QR, vector const &b, vector &x)
C++ version of gsl_linalg_QR_Rsolve().
Definition: linalg.hpp:449
int LQ_decomp(matrix &A, vector &tau)
C++ version of gsl_linalg_LQ_decomp().
Definition: linalg.hpp:646
int bidiag_unpack_B(matrix const &A, vector &diag, vector &superdiag)
C++ version of gsl_linalg_bidiag_unpack_B().
Definition: linalg.hpp:1092
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().
Definition: linalg.hpp:1845
int QRPT_rcond(matrix const &QR, double &rcond, vector &work)
C++ version of gsl_linalg_QRPT_rcond().
Definition: linalg.hpp:1256
int QR_UU_decomp(matrix &U, matrix &S, matrix &T)
C++ version of gsl_linalg_QR_UU_decomp().
Definition: linalg.hpp:2103
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().
Definition: linalg.hpp:788
int QL_unpack(matrix const &QL, vector const &tau, matrix &Q, matrix &L)
C++ version of gsl_linalg_QL_unpack().
Definition: linalg.hpp:2154
int SV_decomp(matrix &A, matrix &V, vector &S, vector &work)
C++ version of gsl_linalg_SV_decomp().
Definition: linalg.hpp:192
int householder_hm1(double tau, matrix &A)
C++ version of gsl_linalg_householder_hm1().
Definition: linalg.hpp:93
int cholesky_band_solve(matrix const &LLT, vector const &b, vector &x)
C++ version of gsl_linalg_cholesky_band_solve().
Definition: linalg.hpp:1595
int cholesky_solve_mat(matrix const &cholesky, matrix const &B, matrix &X)
C++ version of gsl_linalg_cholesky_solve_mat().
Definition: linalg.hpp:1387
int LQ_svx_T(matrix const &LQ, vector const &tau, vector &x)
C++ version of gsl_linalg_LQ_svx_T().
Definition: linalg.hpp:665
int QR_unpack(matrix const &QR, vector const &tau, matrix &Q, matrix &R)
C++ version of gsl_linalg_QR_unpack().
Definition: linalg.hpp:503
int complex_cholesky_solve(matrix_complex const &cholesky, vector_complex const &b, vector_complex &x)
C++ version of gsl_linalg_complex_cholesky_solve().
Definition: linalg.hpp:920
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().
Definition: linalg.hpp:573
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().
Definition: linalg.hpp:1897
int cholesky_decomp_unit(matrix &A, vector &D)
C++ version of gsl_linalg_cholesky_decomp_unit().
Definition: linalg.hpp:904
int LQ_vecQ(matrix const &LQ, vector const &tau, vector &v)
C++ version of gsl_linalg_LQ_vecQ().
Definition: linalg.hpp:712
int QRPT_update(matrix &Q, matrix &R, permutation const &p, vector &u, vector const &v)
C++ version of gsl_linalg_QRPT_update().
Definition: linalg.hpp:638
int complex_tri_invert(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, matrix_complex &T)
C++ version of gsl_linalg_complex_tri_invert().
Definition: linalg.hpp:1800
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().
Definition: linalg.hpp:1048
int complex_LU_invx(matrix_complex &LU, permutation const &p)
C++ version of gsl_linalg_complex_LU_invx().
Definition: linalg.hpp:1172
int LQ_solve_T(matrix const &LQ, vector const &tau, vector const &b, vector &x)
C++ version of gsl_linalg_LQ_solve_T().
Definition: linalg.hpp:656
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().
Definition: linalg.hpp:335
int cholesky_decomp1(matrix &A)
C++ version of gsl_linalg_cholesky_decomp1().
Definition: linalg.hpp:1379
int QR_UR_decomp(matrix &S, matrix &A, matrix &T)
C++ version of gsl_linalg_QR_UR_decomp().
Definition: linalg.hpp:2094
int tri_UL(matrix &LU)
C++ version of gsl_linalg_tri_UL().
Definition: linalg.hpp:1813
int pcholesky_rcond(matrix const &LDLT, permutation const &p, double &rcond, vector &work)
C++ version of gsl_linalg_pcholesky_rcond().
Definition: linalg.hpp:1527
int HH_svx(matrix &A, vector &x)
C++ version of gsl_linalg_HH_svx().
Definition: linalg.hpp:1004
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().
Definition: linalg.hpp:357
int QR_unpack_r(matrix const &QR, matrix const &T, matrix &Q, matrix &R)
C++ version of gsl_linalg_QR_unpack_r().
Definition: linalg.hpp:1201
int ldlt_band_unpack(matrix const &LDLT, matrix &L, vector &D)
C++ version of gsl_linalg_ldlt_band_unpack().
Definition: linalg.hpp:1694
int tri_upper_invert(matrix &T)
C++ version of gsl_linalg_tri_upper_invert().
Definition: linalg.hpp:1762
int LQ_Lsvx_T(matrix const &LQ, vector &x)
C++ version of gsl_linalg_LQ_Lsvx_T().
Definition: linalg.hpp:694
int cholesky_decomp2(matrix &A, vector &S)
C++ version of gsl_linalg_cholesky_decomp2().
Definition: linalg.hpp:1419
int SV_decomp_mod(matrix &A, matrix &X, matrix &V, vector &S, vector &work)
C++ version of gsl_linalg_SV_decomp_mod().
Definition: linalg.hpp:203
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().
Definition: linalg.hpp:2115
int cholesky_svx(matrix const &cholesky, vector &x)
C++ version of gsl_linalg_cholesky_svx().
Definition: linalg.hpp:890
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().
Definition: linalg.hpp:1350
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().
Definition: linalg.hpp:1237
int QR_decomp(matrix &A, vector &tau)
C++ version of gsl_linalg_QR_decomp().
Definition: linalg.hpp:400
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().
Definition: linalg.hpp:1870
int LQ_Lsolve_T(matrix const &LQ, vector const &b, vector &x)
C++ version of gsl_linalg_LQ_Lsolve_T().
Definition: linalg.hpp:686
int QR_QTmat(matrix const &QR, vector const &tau, matrix &A)
C++ version of gsl_linalg_QR_QTmat().
Definition: linalg.hpp:493
int pcholesky_decomp2(matrix &A, permutation &p, vector &S)
C++ version of gsl_linalg_pcholesky_decomp2().
Definition: linalg.hpp:1485
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().
Definition: linalg.hpp:1071
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().
Definition: linalg.hpp:1268
int complex_QR_Qvec(matrix_complex const &QR, vector_complex const &tau, vector_complex &v)
C++ version of gsl_linalg_complex_QR_Qvec().
Definition: linalg.hpp:2013
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().
Definition: linalg.hpp:1337
double householder_transform2(double &alpha, vector &v)
C++ version of gsl_linalg_householder_transform2().
Definition: linalg.hpp:51
int SV_decomp_jacobi(matrix &A, matrix &Q, vector &S)
C++ version of gsl_linalg_SV_decomp_jacobi().
Definition: linalg.hpp:212
int complex_cholesky_svx(matrix_complex const &cholesky, vector_complex &x)
C++ version of gsl_linalg_complex_cholesky_svx().
Definition: linalg.hpp:929
int tri_upper_unit_invert(matrix &T)
C++ version of gsl_linalg_tri_upper_unit_invert().
Definition: linalg.hpp:1775
int tri_LTL(matrix &L)
C++ version of gsl_linalg_tri_LTL().
Definition: linalg.hpp:1807
int symmtd_unpack(matrix const &A, vector const &tau, matrix &Q, vector &diag, vector &subdiag)
C++ version of gsl_linalg_symmtd_unpack().
Definition: linalg.hpp:948
int hermtd_unpack_T(matrix_complex const &A, vector &diag, vector &subdiag)
C++ version of gsl_linalg_hermtd_unpack_T().
Definition: linalg.hpp:987
int QR_QTvec(matrix const &QR, vector const &tau, vector &v)
C++ version of gsl_linalg_QR_QTvec().
Definition: linalg.hpp:475
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().
Definition: linalg.hpp:2060
int hessenberg_unpack_accum(matrix &H, vector &tau, matrix &U)
C++ version of gsl_linalg_hessenberg_unpack_accum().
Definition: linalg.hpp:146
int hermtd_decomp(matrix_complex &A, vector_complex &tau)
C++ version of gsl_linalg_hermtd_decomp().
Definition: linalg.hpp:966
int QL_decomp(matrix &A, vector &tau)
C++ version of gsl_linalg_QL_decomp().
Definition: linalg.hpp:2144
int PTLQ_svx_T(matrix const &LQ, vector const &tau, permutation const &p, vector &x)
C++ version of gsl_linalg_PTLQ_svx_T().
Definition: linalg.hpp:825
complex complex_householder_transform(vector_complex &v)
C++ version of gsl_linalg_complex_householder_transform().
Definition: linalg.hpp:58
int QR_UU_QTvec(matrix const &Y, matrix const &T, vector &b, vector &work)
C++ version of gsl_linalg_QR_UU_QTvec().
Definition: linalg.hpp:2127
int cholesky_band_rcond(matrix const &LLT, double &rcond, vector &work)
C++ version of gsl_linalg_cholesky_band_rcond().
Definition: linalg.hpp:1628
int cholesky_decomp(matrix &A)
C++ version of gsl_linalg_cholesky_decomp().
Definition: linalg.hpp:874
int QR_update(matrix &Q, matrix &R, vector &w, vector const &v)
C++ version of gsl_linalg_QR_update().
Definition: linalg.hpp:466
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().
Definition: linalg.hpp:274
int householder_hv(double tau, vector const &v, vector &w)
C++ version of gsl_linalg_householder_hv().
Definition: linalg.hpp:85
int LQ_LQsolve(matrix &Q, matrix &L, vector const &b, vector &x)
C++ version of gsl_linalg_LQ_LQsolve().
Definition: linalg.hpp:751
double LU_det(matrix &LU, int signum)
C++ version of gsl_linalg_LU_det().
Definition: linalg.hpp:292
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().
Definition: linalg.hpp:1933
int complex_householder_mh(complex &tau, vector_complex const &v, matrix_complex &A)
C++ version of gsl_linalg_complex_householder_mh().
Definition: linalg.hpp:111
int ldlt_rcond(matrix const &LDLT, double &rcond, vector &work)
C++ version of gsl_linalg_ldlt_rcond().
Definition: linalg.hpp:1660
int complex_tri_LHL(matrix_complex &L)
C++ version of gsl_linalg_complex_tri_LHL().
Definition: linalg.hpp:1819
int complex_LU_svx(matrix_complex const &LU, permutation const &p, vector_complex &x)
C++ version of gsl_linalg_complex_LU_svx().
Definition: linalg.hpp:345
int QR_rcond(matrix const &QR, double &rcond, vector &work)
C++ version of gsl_linalg_QR_rcond().
Definition: linalg.hpp:1210
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().
Definition: linalg.hpp:1944
int tri_rcond(CBLAS_UPLO_t Uplo, matrix const &A, double &rcond, vector &work)
C++ version of gsl_linalg_tri_rcond().
Definition: linalg.hpp:1722
int mcholesky_decomp(matrix &A, permutation &p, vector &E)
C++ version of gsl_linalg_mcholesky_decomp().
Definition: linalg.hpp:1537
int pcholesky_invert(matrix const &LDLT, permutation const &p, matrix &Ainv)
C++ version of gsl_linalg_pcholesky_invert().
Definition: linalg.hpp:1517
int LQ_update(matrix &Q, matrix &R, vector const &v, vector &w)
C++ version of gsl_linalg_LQ_update().
Definition: linalg.hpp:741
int symmtd_unpack_T(matrix const &A, vector &diag, vector &subdiag)
C++ version of gsl_linalg_symmtd_unpack_T().
Definition: linalg.hpp:958
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().
Definition: linalg.hpp:836
int ldlt_band_solve(matrix const &LDLT, vector const &b, vector &x)
C++ version of gsl_linalg_ldlt_band_solve().
Definition: linalg.hpp:1677
int QRPT_Rsolve(matrix const &QR, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_QRPT_Rsolve().
Definition: linalg.hpp:618
int mcholesky_invert(matrix const &LDLT, permutation const &p, matrix &Ainv)
C++ version of gsl_linalg_mcholesky_invert().
Definition: linalg.hpp:1577
int SV_solve(matrix const &U, matrix const &Q, vector const &S, vector const &b, vector &x)
C++ version of gsl_linalg_SV_solve().
Definition: linalg.hpp:223
int QRPT_QRsolve(matrix const &Q, matrix const &R, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_QRPT_QRsolve().
Definition: linalg.hpp:607
int symmtd_decomp(matrix &A, vector &tau)
C++ version of gsl_linalg_symmtd_decomp().
Definition: linalg.hpp:937
int LU_invx(matrix &LU, permutation const &p)
C++ version of gsl_linalg_LU_invx().
Definition: linalg.hpp:1164
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().
Definition: linalg.hpp:814
int QR_Rsvx(matrix const &QR, vector &x)
C++ version of gsl_linalg_QR_Rsvx().
Definition: linalg.hpp:457
int QRPT_decomp(matrix &A, vector &tau, permutation &p, int &signum, vector &norm)
C++ version of gsl_linalg_QRPT_decomp().
Definition: linalg.hpp:544
int cholesky_band_svx(matrix const &LLT, vector &x)
C++ version of gsl_linalg_cholesky_band_svx().
Definition: linalg.hpp:1603
int complex_householder_hv(complex &tau, vector_complex const &v, vector_complex &w)
C++ version of gsl_linalg_complex_householder_hv().
Definition: linalg.hpp:120
int R_solve(matrix const &R, vector const &b, vector &x)
C++ version of gsl_linalg_R_solve().
Definition: linalg.hpp:512
int complex_LU_invert(matrix_complex const &LU, permutation const &p, matrix_complex &inverse)
C++ version of gsl_linalg_complex_LU_invert().
Definition: linalg.hpp:368
int hesstri_decomp(matrix &A, matrix &B, matrix &U, matrix &V, vector &work)
C++ version of gsl_linalg_hesstri_decomp().
Definition: linalg.hpp:182
int QR_UD_decomp(matrix &U, vector const &D, matrix &Y, matrix &T)
C++ version of gsl_linalg_QR_UD_decomp().
Definition: linalg.hpp:2071
int QR_matQ(matrix const &QR, vector const &tau, matrix &A)
C++ version of gsl_linalg_QR_matQ().
Definition: linalg.hpp:1191
size_t QRPT_rank(matrix const &QR, double const tol)
C++ version of gsl_linalg_QRPT_rank().
Definition: linalg.hpp:1247
int QR_solve(matrix const &QR, vector const &tau, vector const &b, vector &x)
C++ version of gsl_linalg_QR_solve().
Definition: linalg.hpp:409
int QRPT_solve(matrix const &QR, vector const &tau, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_QRPT_solve().
Definition: linalg.hpp:585
int QR_decomp_old(matrix &A, vector &tau)
C++ version of gsl_linalg_QR_decomp_old().
Definition: linalg.hpp:1907
int cholesky_invert(matrix &cholesky)
C++ version of gsl_linalg_cholesky_invert().
Definition: linalg.hpp:897
int HH_solve(matrix &A, vector const &b, vector &x)
C++ version of gsl_linalg_HH_solve().
Definition: linalg.hpp:996
int L_solve_T(matrix const &L, vector const &b, vector &x)
C++ version of gsl_linalg_L_solve_T().
Definition: linalg.hpp:703
int complex_householder_hm(complex &tau, vector_complex const &v, matrix_complex &A)
C++ version of gsl_linalg_complex_householder_hm().
Definition: linalg.hpp:102
int hessenberg_unpack(matrix &H, vector &tau, matrix &U)
C++ version of gsl_linalg_hessenberg_unpack().
Definition: linalg.hpp:137
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().
Definition: linalg.hpp:1036
int pcholesky_svx(matrix const &LDLT, permutation const &p, vector &x)
C++ version of gsl_linalg_pcholesky_svx().
Definition: linalg.hpp:1476
int mcholesky_solve(matrix const &LDLT, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_mcholesky_solve().
Definition: linalg.hpp:1547
int cholesky_rcond(matrix const &LLT, double &rcond, vector &work)
C++ version of gsl_linalg_cholesky_rcond().
Definition: linalg.hpp:1448
int SV_leverage(matrix const &U, vector &h)
C++ version of gsl_linalg_SV_leverage().
Definition: linalg.hpp:1124
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().
Definition: linalg.hpp:1966
int LQ_unpack(matrix const &LQ, vector const &tau, matrix &Q, matrix &L)
C++ version of gsl_linalg_LQ_unpack().
Definition: linalg.hpp:731
int QR_UZ_decomp(matrix &S, matrix &A, matrix &T)
C++ version of gsl_linalg_QR_UZ_decomp().
Definition: linalg.hpp:2136
int PTLQ_decomp(matrix &A, vector &tau, permutation &p, int &signum, vector &norm)
C++ version of gsl_linalg_PTLQ_decomp().
Definition: linalg.hpp:775
int tri_lower_rcond(matrix const &A, double &rcond, vector &work)
C++ version of gsl_linalg_tri_lower_rcond().
Definition: linalg.hpp:1740
int LQ_vecQT(matrix const &LQ, vector const &tau, vector &v)
C++ version of gsl_linalg_LQ_vecQT().
Definition: linalg.hpp:721
int complex_QR_svx(matrix_complex const &QR, vector_complex const &tau, vector_complex &x)
C++ version of gsl_linalg_complex_QR_svx().
Definition: linalg.hpp:1954
int PTLQ_Lsolve_T(matrix const &LQ, permutation const &p, vector const &b, vector &x)
C++ version of gsl_linalg_PTLQ_Lsolve_T().
Definition: linalg.hpp:847
int LU_svx(matrix const &LU, permutation const &p, vector &x)
C++ version of gsl_linalg_LU_svx().
Definition: linalg.hpp:262
int PTLQ_update(matrix &Q, matrix &L, permutation const &p, vector const &v, vector &w)
C++ version of gsl_linalg_PTLQ_update().
Definition: linalg.hpp:867
int QR_lssolve(matrix const &QR, vector const &tau, vector const &b, vector &x, vector &residual)
C++ version of gsl_linalg_QR_lssolve().
Definition: linalg.hpp:429
double householder_transform(vector &v)
C++ version of gsl_linalg_householder_transform().
Definition: linalg.hpp:43
int ldlt_decomp(matrix &A)
C++ version of gsl_linalg_ldlt_decomp().
Definition: linalg.hpp:1635
int solve_tridiag(vector const &diag, vector const &abovediag, vector const &belowdiag, vector const &b, vector &x)
C++ version of gsl_linalg_solve_tridiag().
Definition: linalg.hpp:1025
int PTLQ_Lsvx_T(matrix const &LQ, permutation const &p, vector &x)
C++ version of gsl_linalg_PTLQ_Lsvx_T().
Definition: linalg.hpp:856
int hessenberg_decomp(matrix &A, vector &tau)
C++ version of gsl_linalg_hessenberg_decomp().
Definition: linalg.hpp:128
int complex_cholesky_decomp(matrix_complex &A)
C++ version of gsl_linalg_complex_cholesky_decomp().
Definition: linalg.hpp:911
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().
Definition: linalg.hpp:2083
int complex_tri_UL(matrix_complex &LU)
C++ version of gsl_linalg_complex_tri_UL().
Definition: linalg.hpp:1826
int pcholesky_svx2(matrix const &LDLT, permutation const &p, vector const &S, vector &x)
C++ version of gsl_linalg_pcholesky_svx2().
Definition: linalg.hpp:1507
int LU_invert(matrix const &LU, permutation const &p, matrix &inverse)
C++ version of gsl_linalg_LU_invert().
Definition: linalg.hpp:284
int bidiag_unpack2(matrix &A, vector &tau_U, vector &tau_V, matrix &V)
C++ version of gsl_linalg_bidiag_unpack2().
Definition: linalg.hpp:1083
int householder_mh(double tau, vector const &v, matrix &A)
C++ version of gsl_linalg_householder_mh().
Definition: linalg.hpp:76
int householder_left(double const tau, vector const &v, matrix &A, vector &work)
C++ version of gsl_linalg_householder_left().
Definition: linalg.hpp:1134
double complex_LU_lndet(matrix_complex &LU)
C++ version of gsl_linalg_complex_LU_lndet().
Definition: linalg.hpp:384
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().
Definition: linalg.hpp:1857
int householder_right(double const tau, vector const &v, matrix &A, vector &work)
C++ version of gsl_linalg_householder_right().
Definition: linalg.hpp:1145
int QR_QTmat_r(matrix const &QR, matrix const &T, matrix &B, matrix &work)
C++ version of gsl_linalg_QR_QTmat_r().
Definition: linalg.hpp:1182
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().
Definition: linalg.hpp:2047
int hessenberg_submatrix(matrix &M, matrix &A, size_t top, vector &tau)
C++ version of gsl_linalg_hessenberg_submatrix().
Definition: linalg.hpp:162
int complex_LU_decomp(matrix_complex &A, permutation &p, int &signum)
C++ version of gsl_linalg_complex_LU_decomp().
Definition: linalg.hpp:325
int QR_svx(matrix const &QR, vector const &tau, vector &x)
C++ version of gsl_linalg_QR_svx().
Definition: linalg.hpp:418
int householder_hm(double tau, vector const &v, matrix &A)
C++ version of gsl_linalg_householder_hm().
Definition: linalg.hpp:67
int cholesky_solve2(matrix const &LLT, vector const &S, vector const &b, vector &x)
C++ version of gsl_linalg_cholesky_solve2().
Definition: linalg.hpp:1438
int QR_QRsolve(matrix &Q, matrix &R, vector const &b, vector &x)
C++ version of gsl_linalg_QR_QRsolve().
Definition: linalg.hpp:440
int cholesky_band_decomp(matrix &A)
C++ version of gsl_linalg_cholesky_band_decomp().
Definition: linalg.hpp:1586
int complex_QR_decomp_r(matrix_complex &A, matrix_complex &T)
C++ version of gsl_linalg_complex_QR_decomp_r().
Definition: linalg.hpp:1923
int ldlt_solve(matrix const &LDLT, vector const &b, vector &x)
C++ version of gsl_linalg_ldlt_solve().
Definition: linalg.hpp:1643
int LQ_lssolve(matrix const &LQ, vector const &tau, vector const &b, vector &x, vector &residual)
C++ version of gsl_linalg_LQ_lssolve().
Definition: linalg.hpp:1362
void givens(double const a, double const b, double &c, double &s)
C++ version of gsl_linalg_givens().
Definition: linalg.hpp:1835
int QRPT_Rsvx(matrix const &QR, permutation const &p, vector &x)
C++ version of gsl_linalg_QRPT_Rsvx().
Definition: linalg.hpp:627
int LU_sgndet(matrix &lu, int signum)
C++ version of gsl_linalg_LU_sgndet().
Definition: linalg.hpp:306
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().
Definition: linalg.hpp:1319
int ldlt_band_svx(matrix const &LDLT, vector &x)
C++ version of gsl_linalg_ldlt_band_svx().
Definition: linalg.hpp:1685
int balance_accum(matrix &A, vector &D)
C++ version of gsl_linalg_balance_accum().
Definition: linalg.hpp:1108
int cholesky_scale(matrix const &A, vector &S)
C++ version of gsl_linalg_cholesky_scale().
Definition: linalg.hpp:1403
int tri_upper_rcond(matrix const &A, double &rcond, vector &work)
C++ version of gsl_linalg_tri_upper_rcond().
Definition: linalg.hpp:1731
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().
Definition: linalg.hpp:676
int tri_lower_unit_invert(matrix &T)
C++ version of gsl_linalg_tri_lower_unit_invert().
Definition: linalg.hpp:1782
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().
Definition: linalg.hpp:977
int LQ_QTvec(matrix const &LQ, vector const &tau, vector &v)
C++ version of gsl_linalg_LQ_QTvec().
Definition: linalg.hpp:1372
int pcholesky_solve2(matrix const &LDLT, permutation const &p, vector const &S, vector const &b, vector &x)
C++ version of gsl_linalg_pcholesky_solve2().
Definition: linalg.hpp:1496
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().
Definition: linalg.hpp:1155
int tri_lower_invert(matrix &T)
C++ version of gsl_linalg_tri_lower_invert().
Definition: linalg.hpp:1769
int ldlt_band_decomp(matrix &A)
C++ version of gsl_linalg_ldlt_band_decomp().
Definition: linalg.hpp:1668
double const E
e
Definition: math.hpp:36
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().
Definition: sf_ellint.hpp:203
double b(int order, double qq)
C++ version of gsl_sf_mathieu_b().
Definition: sf_mathieu.hpp:298
double a(int order, double qq)
C++ version of gsl_sf_mathieu_a().
Definition: sf_mathieu.hpp:272
The gsl package creates an interface to the GNU Scientific Library for C++.
Definition: blas.hpp:34