- Timestamp:
- Feb 17, 2007, 1:25:36 PM (16 years ago)
- Location:
- trunk/yat
- Files:
-
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/yat/regression/MultiDimensional.h
r739 r751 56 56 const utility::matrix& covariance(void) const; 57 57 58 /// 59 /// Function fitting parameters of the linear model by miminizing 60 /// the quadratic deviation between model and data. 61 /// 58 /** 59 \brief Function fitting parameters of the linear model by 60 miminizing the quadratic deviation between model and data. 61 62 \throw A GSL_error exception is thrown if memory allocation 63 fails or the underlying GSL calls fails (usually matrix 64 dimension errors). 65 */ 62 66 void fit(const utility::matrix& X, const utility::vector& y); 63 67 -
trunk/yat/regression/MultiDimensionalWeighted.h
r741 r751 56 56 double chisq(void) const; 57 57 58 /// 59 /// @see gsl_multifit_wlinear 60 /// 58 /** 59 \see gsl_multifit_wlinear 60 61 \throw A GSL_error exception is thrown if memory allocation 62 fails or the underlying GSL calls fails (usually matrix 63 dimension errors). 64 */ 61 65 void fit(const utility::matrix& X, const utility::vector& y, 62 66 const utility::vector& w); -
trunk/yat/utility/Exception.h
r750 r751 53 53 54 54 inline GSL_error(std::string message, int gsl_status) throw() 55 : std::runtime_error("GSL_error: " + message + gsl_strerror(gsl_status)) {} 55 : std::runtime_error("GSL_error: " + message + " " + 56 gsl_strerror(gsl_status)) {} 56 57 }; 57 58 -
trunk/yat/utility/SVD.cc
r719 r751 23 23 24 24 #include "SVD.h" 25 #include <sstream> 25 26 26 27 namespace theplu { … … 40 41 41 42 42 intSVD::decompose(SVDalgorithm algo)43 void SVD::decompose(SVDalgorithm algo) 43 44 { 45 int status=0; 44 46 switch (algo) { 45 case GolubReinsch: 46 return golub_reinsch(); 47 case ModifiedGolubReinsch: 48 return modified_golub_reinsch(); 49 case Jacobi: 50 return jacobi(); 47 case GolubReinsch: 48 status=golub_reinsch(); 49 break; 50 case ModifiedGolubReinsch: 51 status=modified_golub_reinsch(); 52 break; 53 case Jacobi: 54 status=jacobi(); 55 break; 56 default: 57 std::stringstream ss("SVD::decompose "); 58 ss << algo << " is not a valid SVDalgorithm"; 59 throw GSL_error(ss.str()); 51 60 } 52 // Jari, the program should never end up here, return values should be 53 // something different from normal GSL return values, or maybe 54 // throw an exception. 55 return 0; 61 if (status) 62 throw utility::GSL_error(std::string("SVD::decompose",status)); 56 63 } 57 64 … … 88 95 89 96 90 intSVD::solve(const utility::vector& b, utility::vector& x)97 void SVD::solve(const utility::vector& b, utility::vector& x) 91 98 { 92 return gsl_linalg_SV_solve(U_.gsl_matrix_p(), V_.gsl_matrix_p(), 93 s_.gsl_vector_p(), b.gsl_vector_p(), 94 x.gsl_vector_p()); 99 int status=gsl_linalg_SV_solve(U_.gsl_matrix_p(), V_.gsl_matrix_p(), 100 s_.gsl_vector_p(), b.gsl_vector_p(), 101 x.gsl_vector_p()); 102 if (status) 103 throw utility::GSL_error(std::string("SVD::solve",status)); 95 104 } 96 105 -
trunk/yat/utility/SVD.h
r719 r751 34 34 namespace utility { 35 35 36 // Jari check that interface is complete37 38 36 /** 39 37 Class encapsulating GSL methods for singular value decomposition, … … 47 45 V = Orthogonal matrix, size NxN\n 48 46 */ 49 50 47 class SVD 51 48 { 52 49 public: 53 50 54 / //55 ///A number of SVD algorithms are implemented in GSL. They have56 ///their strengths and weaknesses, check the GSL documentation.57 ///58 ///There are restrictions on the matrix dimensions depending on59 /// which SVD algorithm is used. From the GSL's SVD source code60 /// one finds that the Golub-Reinsch algorithm implementation will61 /// not work on matrices with fewer rows than columns, the same is62 /// alsotrue for the modified Golub-Reinsch algorithm.63 ///64 /// @see GSL's SVD documentation.65 ///51 /** 52 A number of SVD algorithms are implemented in GSL. They have 53 their strengths and weaknesses, check the GSL documentation. 54 55 There are restrictions on the matrix dimensions depending on 56 which SVD algorithm is used. From the GSL's SVD source code one 57 finds that the Golub-Reinsch algorithm implementation will not 58 work on matrices with fewer rows than columns, the same is also 59 true for the modified Golub-Reinsch algorithm. 60 61 \see GSL's SVD documentation. 62 */ 66 63 enum SVDalgorithm { 67 64 GolubReinsch, … … 70 67 }; 71 68 72 /// 73 /// Constructs an SVD object using the matrix Ain as only input. The 74 /// input matrix is copied for further use in the object. 75 /// 69 /** 70 \brief Constructs an SVD object using the matrix Ain as only 71 input. The input matrix is copied for further use in the 72 object. 73 */ 76 74 SVD(const utility::matrix& Ain); 77 75 78 / //79 /// @brief The destructor80 ///76 /** 77 \brief The destructor 78 */ 81 79 ~SVD(void); 82 80 83 /// 84 /// This function will perform SVD with the method specified by \a 85 /// algo. 86 /// 87 /// @return Whatever GSL returns. 88 /// 89 int decompose(SVDalgorithm algo=GolubReinsch); 81 /** 82 \brief This function will perform SVD with the method specified 83 by \a algo. 90 84 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 /// 85 \throw GSL_error if the underlying GSL function fails. 86 */ 87 void decompose(SVDalgorithm algo=GolubReinsch); 88 89 /** 90 \brief Access to the s vector. 91 92 \return A copy of the s vector. 93 94 \note If decompose() has not been run the outcome of the call 95 is undefined. 96 */ 99 97 const utility::vector& s(void) const; 100 98 101 /// 102 /// @brief Solve the system \f$ Ax=b \f$ using the decomposition 103 /// of A. 104 /// 105 /// @note If decompose() has not been run the outcome of the call 106 /// is undefined. 107 /// 108 /// @return Whatever GSL returns. 109 /// 110 int solve(const utility::vector& b, utility::vector& x); 99 /** 100 \brief Solve the system \f$ Ax=b \f$ using the decomposition of 101 A. 111 102 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 /// 103 \note If decompose() has not been run the outcome of the call 104 is undefined. 105 106 \throw GSL_error if the underlying GSL function fails. 107 */ 108 void solve(const utility::vector& b, utility::vector& x); 109 110 /** 111 \brief Access to the U matrix. 112 113 \return A copy of the U matrix. 114 115 \note If decompose() has not been run the outcome of the call 116 is undefined. 117 */ 120 118 const utility::matrix& U(void) const; 121 119 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 call128 ///is undefined.129 ///120 /** 121 \brief Access to the V matrix. 122 123 \return A copy of the V matrix. 124 125 \note If decompose() has not been run the outcome of the call 126 is undefined. 127 */ 130 128 const utility::matrix& V(void) const; 131 129 132 130 private: 131 /** 132 \brief Call GSL's Jacobi algorithm for SVD. 133 134 \return gsl_error status. 135 */ 133 136 int jacobi(void); 137 138 /** 139 \brief Call GSL's Golub-Reinsch algorithm for SVD. 140 141 \return gsl_error status. 142 */ 134 143 int golub_reinsch(void); 144 145 /** 146 \brief Call GSL's modified Golub-Reinch algorithm for SVD. 147 148 \return gsl_error status. 149 */ 135 150 int modified_golub_reinsch(void); 136 151
Note: See TracChangeset
for help on using the changeset viewer.