1 | // $Id: MultiDimensionalWeighted.cc 1392 2008-07-28 19:35:30Z peter $ |
---|
2 | |
---|
3 | /* |
---|
4 | Copyright (C) 2006, 2007, 2008 Jari Häkkinen, Peter Johansson |
---|
5 | |
---|
6 | This file is part of the yat library, http://dev.thep.lu.se/yat |
---|
7 | |
---|
8 | The yat library is free software; you can redistribute it and/or |
---|
9 | modify it under the terms of the GNU General Public License as |
---|
10 | published by the Free Software Foundation; either version 2 of the |
---|
11 | License, or (at your option) any later version. |
---|
12 | |
---|
13 | The yat library is distributed in the hope that it will be useful, |
---|
14 | but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
15 | MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
---|
16 | General Public License for more details. |
---|
17 | |
---|
18 | You should have received a copy of the GNU General Public License |
---|
19 | along with this program; if not, write to the Free Software |
---|
20 | Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA |
---|
21 | 02111-1307, USA. |
---|
22 | */ |
---|
23 | |
---|
24 | #include "MultiDimensionalWeighted.h" |
---|
25 | #include "yat/statistics/AveragerWeighted.h" |
---|
26 | #include "yat/utility/Matrix.h" |
---|
27 | #include "yat/utility/Vector.h" |
---|
28 | |
---|
29 | #include <cassert> |
---|
30 | |
---|
31 | namespace theplu { |
---|
32 | namespace yat { |
---|
33 | namespace regression { |
---|
34 | |
---|
35 | MultiDimensionalWeighted::MultiDimensionalWeighted(void) |
---|
36 | : chisquare_(0), work_(NULL) |
---|
37 | { |
---|
38 | } |
---|
39 | |
---|
40 | MultiDimensionalWeighted::~MultiDimensionalWeighted(void) |
---|
41 | { |
---|
42 | if (work_) |
---|
43 | gsl_multifit_linear_free(work_); |
---|
44 | } |
---|
45 | |
---|
46 | |
---|
47 | double MultiDimensionalWeighted::chisq() const |
---|
48 | { |
---|
49 | return chisquare_; |
---|
50 | } |
---|
51 | |
---|
52 | |
---|
53 | void MultiDimensionalWeighted::fit(const utility::Matrix& x, |
---|
54 | const utility::VectorBase& y, |
---|
55 | const utility::VectorBase& w) |
---|
56 | { |
---|
57 | assert(y.size()==w.size()); |
---|
58 | assert(x.rows()==y.size()); |
---|
59 | |
---|
60 | covariance_.resize(x.columns(),x.columns()); |
---|
61 | fit_parameters_ = utility::Vector(x.columns()); |
---|
62 | if (work_) |
---|
63 | gsl_multifit_linear_free(work_); |
---|
64 | if (!(work_=gsl_multifit_linear_alloc(x.rows(),fit_parameters_.size()))) |
---|
65 | throw utility::GSL_error("MultiDimensionalWeighted::fit failed to allocate memory"); |
---|
66 | int status = gsl_multifit_wlinear(x.gsl_matrix_p(), w.gsl_vector_p(), |
---|
67 | y.gsl_vector_p(), |
---|
68 | fit_parameters_.gsl_vector_p(), |
---|
69 | covariance_.gsl_matrix_p(), &chisquare_, |
---|
70 | work_); |
---|
71 | if (status) |
---|
72 | throw utility::GSL_error(std::string("MultiDimensionalWeighted::fit", |
---|
73 | status)); |
---|
74 | |
---|
75 | statistics::AveragerWeighted aw; |
---|
76 | add(aw, y.begin(), y.end(), w.begin()); |
---|
77 | s2_ = chisquare_ / (aw.n()-fit_parameters_.size()); |
---|
78 | covariance_ *= s2_; |
---|
79 | } |
---|
80 | |
---|
81 | |
---|
82 | const utility::Vector& MultiDimensionalWeighted::fit_parameters(void) const |
---|
83 | { |
---|
84 | return fit_parameters_; |
---|
85 | } |
---|
86 | |
---|
87 | |
---|
88 | double MultiDimensionalWeighted::predict(const utility::VectorBase& x) const |
---|
89 | { |
---|
90 | assert(x.size()==fit_parameters_.size()); |
---|
91 | return fit_parameters_ * x; |
---|
92 | } |
---|
93 | |
---|
94 | |
---|
95 | double MultiDimensionalWeighted::prediction_error2(const utility::VectorBase& x, |
---|
96 | const double w) const |
---|
97 | { |
---|
98 | return standard_error2(x) + s2(w); |
---|
99 | } |
---|
100 | |
---|
101 | |
---|
102 | double MultiDimensionalWeighted::s2(const double w) const |
---|
103 | { |
---|
104 | return s2_/w; |
---|
105 | } |
---|
106 | |
---|
107 | |
---|
108 | double |
---|
109 | MultiDimensionalWeighted::standard_error2(const utility::VectorBase& x) const |
---|
110 | { |
---|
111 | double c = 0; |
---|
112 | for (size_t i=0; i<x.size(); ++i){ |
---|
113 | c += covariance_(i,i)*x(i)*x(i); |
---|
114 | for (size_t j=i+1; j<x.size(); ++j) |
---|
115 | c += 2*covariance_(i,j)*x(i)*x(j); |
---|
116 | } |
---|
117 | return c; |
---|
118 | } |
---|
119 | |
---|
120 | }}} // of namespaces regression, yat, and theplu |
---|