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

Last change on this file since 767 was 767, checked in by Peter, 15 years ago

Fixes #65

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