source: trunk/yat/utility/SVD.h @ 1486

Last change on this file since 1486 was 1486, checked in by Jari Häkkinen, 13 years ago

Addresses #436.

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 3.7 KB
Line 
1#ifndef _theplu_yat_utility_svd_
2#define _theplu_yat_utility_svd_
3
4// $Id: SVD.h 1486 2008-09-09 21:17:19Z jari $
5
6/*
7  Copyright (C) 2003 Daniel Dalevi, Jari Häkkinen
8  Copyright (C) 2004 Jari Häkkinen
9  Copyright (C) 2005 Jari Häkkinen, Peter Johansson
10  Copyright (C) 2006 Jari Häkkinen, Markus Ringnér
11  Copyright (C) 2007 Jari Häkkinen, Peter Johansson
12  Copyright (C) 2008 Peter Johansson
13
14  This file is part of the yat library, http://dev.thep.lu.se/yat
15
16  The yat library is free software; you can redistribute it and/or
17  modify it under the terms of the GNU General Public License as
18  published by the Free Software Foundation; either version 3 of the
19  License, or (at your option) any later version.
20
21  The yat library is distributed in the hope that it will be useful,
22  but WITHOUT ANY WARRANTY; without even the implied warranty of
23  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
24  General Public License for more details.
25
26  You should have received a copy of the GNU General Public License
27  along with this program; if not, write to the Free Software
28  Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA
29  02111-1307, USA.
30*/
31
32#include "Matrix.h"
33#include "Vector.h"
34
35#include <gsl/gsl_linalg.h>
36
37namespace theplu {
38namespace yat {
39namespace utility {
40
41  class VectorBase;
42
43  /**
44     @brief Singular Value Decomposition
45
46     Class encapsulating GSL methods for singular value
47     decomposition, SVD.
48
49     A = U S V' = (MxN)(NxN)(NxN) = (MxN)\n   
50
51     A = Matrix to be decomposed, size MxN\n
52     U = Orthogonal matrix, size MxN\n
53     S = Diagonal matrix of singular values, size NxN\n
54     V = Orthogonal matrix, size NxN\n
55  */
56  class SVD
57  {
58  public:
59
60    /**
61       A number of SVD algorithms are implemented in GSL. They have
62       their strengths and weaknesses.
63   
64       \see GSL's SVD documentation.
65    */
66    enum SVDalgorithm {
67      GolubReinsch,
68      ModifiedGolubReinsch,
69      Jacobi
70    };
71
72    /**
73       \brief Constructs an SVD object using the matrix Ain as only
74       input. The input matrix is copied for further use in the
75       object.
76
77       \note Number of rows must be equal or larger than number of columns.
78    */
79    SVD(const utility::Matrix& Ain);
80
81    /**
82       \brief The destructor
83    */
84    ~SVD(void);
85
86    /**
87       \brief This function will perform SVD with the method specified
88       by \a algo.
89
90       \throw GSL_error if the underlying GSL function fails.
91    */
92    void decompose(SVDalgorithm algo=GolubReinsch);
93
94    /**
95       \brief Access to the s vector.
96
97       \return A copy of the s vector.
98
99       \note If decompose() has not been run the outcome of the call
100       is undefined.
101    */
102    const utility::Vector& s(void) const;
103
104    /**
105       \brief Solve the system \f$ Ax=b \f$ using the decomposition of
106       A.
107
108       \note If decompose() has not been run the outcome of the call
109       is undefined.
110
111       \throw GSL_error if the underlying GSL function fails.
112    */
113    void solve(const utility::VectorBase& b, utility::Vector& x);
114
115    /**
116       \brief Access to the U matrix.
117
118       \return A copy of the U matrix.
119
120       \note If decompose() has not been run the outcome of the call
121       is undefined.
122    */
123    const utility::Matrix& U(void) const;
124
125    /**
126       \brief Access to the V matrix.
127
128       \return A copy of the V matrix.
129
130       \note If decompose() has not been run the outcome of the call
131       is undefined.
132    */
133    const utility::Matrix& V(void) const;
134
135  private:
136    /**
137       \brief Call GSL's Jacobi algorithm for SVD.
138
139       \return gsl_error status.
140    */
141    int jacobi(void);
142
143    /**
144       \brief Call GSL's Golub-Reinsch algorithm for SVD.
145
146       \return gsl_error status.
147    */
148    int golub_reinsch(void);
149
150    /**
151       \brief Call GSL's modified Golub-Reinch algorithm for SVD.
152
153       \return gsl_error status.
154    */
155    int modified_golub_reinsch(void);
156
157    utility::Matrix U_, V_;
158    utility::Vector s_;
159  }; 
160
161}}} // of namespace utility, yat, and theplu
162
163#endif
Note: See TracBrowser for help on using the repository browser.