1 | // $Id: SVD.cc 3736 2018-05-22 00:59:29Z peter $ |
---|
2 | |
---|
3 | /* |
---|
4 | Copyright (C) 2003 Daniel Dalevi |
---|
5 | Copyright (C) 2004 Jari Häkkinen |
---|
6 | Copyright (C) 2005 Jari Häkkinen, Peter Johansson |
---|
7 | Copyright (C) 2006 Jari Häkkinen |
---|
8 | Copyright (C) 2007, 2008 Jari Häkkinen, Peter Johansson |
---|
9 | Copyright (C) 2009 Jari Häkkinen |
---|
10 | Copyright (C) 2011, 2012 Peter Johansson |
---|
11 | |
---|
12 | This file is part of the yat library, http://dev.thep.lu.se/yat |
---|
13 | |
---|
14 | The yat library is free software; you can redistribute it and/or |
---|
15 | modify it under the terms of the GNU General Public License as |
---|
16 | published by the Free Software Foundation; either version 3 of the |
---|
17 | License, or (at your option) any later version. |
---|
18 | |
---|
19 | The yat library is distributed in the hope that it will be useful, |
---|
20 | but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
21 | MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
---|
22 | General Public License for more details. |
---|
23 | |
---|
24 | You should have received a copy of the GNU General Public License |
---|
25 | along with yat. If not, see <http://www.gnu.org/licenses/>. |
---|
26 | */ |
---|
27 | |
---|
28 | #include <config.h> |
---|
29 | |
---|
30 | #include "SVD.h" |
---|
31 | |
---|
32 | #include "Exception.h" |
---|
33 | #include "Vector.h" |
---|
34 | #include "VectorBase.h" |
---|
35 | |
---|
36 | #include <cassert> |
---|
37 | #include <sstream> |
---|
38 | #include <utility> |
---|
39 | |
---|
40 | namespace theplu { |
---|
41 | namespace yat { |
---|
42 | namespace utility { |
---|
43 | |
---|
44 | |
---|
45 | SVD::SVD(const Matrix& Ain) |
---|
46 | : U_(Ain), V_(Ain.columns(),Ain.columns()), s_(Ain.columns()) |
---|
47 | { |
---|
48 | } |
---|
49 | |
---|
50 | |
---|
51 | #if YAT_HAVE_RVALE |
---|
52 | SVD::SVD(Matrix&& Ain) |
---|
53 | : U_(std::move(Ain)), V_(Ain.columns(),Ain.columns()), s_(Ain.columns()) |
---|
54 | { |
---|
55 | } |
---|
56 | #endif |
---|
57 | |
---|
58 | |
---|
59 | SVD::~SVD(void) |
---|
60 | { |
---|
61 | } |
---|
62 | |
---|
63 | |
---|
64 | void SVD::decompose(SVDalgorithm algo) |
---|
65 | { |
---|
66 | assert(U_.rows()>=U_.columns()); |
---|
67 | int status=0; |
---|
68 | switch (algo) { |
---|
69 | case GolubReinsch: |
---|
70 | status=golub_reinsch(); |
---|
71 | break; |
---|
72 | case ModifiedGolubReinsch: |
---|
73 | status=modified_golub_reinsch(); |
---|
74 | break; |
---|
75 | case Jacobi: |
---|
76 | status=jacobi(); |
---|
77 | break; |
---|
78 | default: |
---|
79 | std::stringstream ss("SVD::decompose "); |
---|
80 | ss << algo << " is not a valid SVDalgorithm"; |
---|
81 | throw GSL_error(ss.str()); |
---|
82 | } |
---|
83 | if (status) |
---|
84 | throw GSL_error(std::string("SVD::decompose",status)); |
---|
85 | } |
---|
86 | |
---|
87 | |
---|
88 | int SVD::golub_reinsch(void) |
---|
89 | { |
---|
90 | Vector w(U_.columns()); |
---|
91 | return gsl_linalg_SV_decomp(U_.gsl_matrix_p(), V_.gsl_matrix_p(), |
---|
92 | s_.gsl_vector_p(), w.gsl_vector_p()); |
---|
93 | } |
---|
94 | |
---|
95 | |
---|
96 | int SVD::jacobi(void) |
---|
97 | { |
---|
98 | return gsl_linalg_SV_decomp_jacobi(U_.gsl_matrix_p(), V_.gsl_matrix_p(), |
---|
99 | s_.gsl_vector_p()); |
---|
100 | } |
---|
101 | |
---|
102 | |
---|
103 | int SVD::modified_golub_reinsch(void) |
---|
104 | { |
---|
105 | Vector w(U_.columns()); |
---|
106 | Matrix X(U_.columns(),U_.columns()); |
---|
107 | return gsl_linalg_SV_decomp_mod(U_.gsl_matrix_p(), X.gsl_matrix_p(), |
---|
108 | V_.gsl_matrix_p(), s_.gsl_vector_p(), |
---|
109 | w.gsl_vector_p()); |
---|
110 | } |
---|
111 | |
---|
112 | |
---|
113 | const Vector& SVD::s(void) const |
---|
114 | { |
---|
115 | return s_; |
---|
116 | } |
---|
117 | |
---|
118 | |
---|
119 | void SVD::solve(const VectorBase& b, Vector& x) |
---|
120 | { |
---|
121 | int status=gsl_linalg_SV_solve(U_.gsl_matrix_p(), V_.gsl_matrix_p(), |
---|
122 | s_.gsl_vector_p(), b.gsl_vector_p(), |
---|
123 | x.gsl_vector_p()); |
---|
124 | if (status) |
---|
125 | throw GSL_error(std::string("SVD::solve",status)); |
---|
126 | } |
---|
127 | |
---|
128 | |
---|
129 | const Matrix& SVD::U(void) const |
---|
130 | { |
---|
131 | return U_; |
---|
132 | } |
---|
133 | |
---|
134 | |
---|
135 | const Matrix& SVD::V(void) const |
---|
136 | { |
---|
137 | return V_; |
---|
138 | } |
---|
139 | |
---|
140 | }}} // of namespace utility, yat, and theplu |
---|