source: trunk/yat/regression/MultiDimensional.cc @ 1615

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

Addresses #436. GPL license copy reference should also be updated.

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 2.8 KB
Line 
1// $Id: MultiDimensional.cc 1487 2008-09-10 08:41:36Z jari $
2
3/*
4  Copyright (C) 2005 Jari Häkkinen
5  Copyright (C) 2006, 2007, 2008 Jari Häkkinen, Peter Johansson
6
7  This file is part of the yat library, http://dev.thep.lu.se/yat
8
9  The yat library is free software; you can redistribute it and/or
10  modify it under the terms of the GNU General Public License as
11  published by the Free Software Foundation; either version 3 of the
12  License, or (at your option) any later version.
13
14  The yat library is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17  General Public License for more details.
18
19  You should have received a copy of the GNU General Public License
20  along with yat. If not, see <http://www.gnu.org/licenses/>.
21*/
22
23#include "MultiDimensional.h"
24#include "yat/utility/Exception.h"
25#include "yat/utility/Matrix.h"
26#include "yat/utility/VectorBase.h"
27#include "yat/utility/Vector.h"
28
29#include <cassert>
30
31namespace theplu {
32namespace yat {
33namespace regression {
34
35
36  MultiDimensional::MultiDimensional(void)
37    : chisquare_(0), work_(NULL)
38  {
39  }
40
41
42  MultiDimensional::~MultiDimensional(void)
43  {
44    if (work_)
45      gsl_multifit_linear_free(work_);
46  }
47
48
49  const utility::Matrix& MultiDimensional::covariance(void) const
50  {
51    return covariance_;
52  }
53
54
55  void MultiDimensional::fit(const utility::Matrix& x, 
56                             const utility::VectorBase& y)
57  {
58    assert(x.rows()==y.size());
59    covariance_.resize(x.columns(),x.columns());
60    fit_parameters_ = utility::Vector(x.columns());
61    if (work_)
62      gsl_multifit_linear_free(work_);
63    if (!(work_=gsl_multifit_linear_alloc(x.rows(),fit_parameters_.size())))
64      throw utility::GSL_error("MultiDimensional::fit failed to allocate memory");
65
66    int status = gsl_multifit_linear(x.gsl_matrix_p(), y.gsl_vector_p(),
67                                     fit_parameters_.gsl_vector_p(),
68                                     covariance_.gsl_matrix_p(), &chisquare_,
69                                     work_);
70    if (status)
71      throw utility::GSL_error(std::string("MultiDimensional::fit",status));
72
73    s2_ = chisquare_/(x.rows()-x.columns());
74  }
75
76  const utility::Vector& MultiDimensional::fit_parameters(void) const
77  { 
78    return fit_parameters_; 
79  }
80
81
82  double MultiDimensional::chisq(void) const
83  {
84    return chisquare_;
85  }
86
87
88  double MultiDimensional::predict(const utility::VectorBase& x) const 
89  {
90    assert(x.size()==fit_parameters_.size());
91    return fit_parameters_ * x;
92  }
93
94
95  double MultiDimensional::prediction_error2(const utility::VectorBase& x) const
96  {
97    return standard_error2(x) + s2_;
98  }
99
100
101  double MultiDimensional::standard_error2(const utility::VectorBase& x) const
102  {
103    double s2 = 0;
104    for (size_t i=0; i<x.size(); ++i){
105      s2 += covariance_(i,i)*x(i)*x(i);
106      for (size_t j=i+1; j<x.size(); ++j)
107        s2 += 2*covariance_(i,j)*x(i)*x(j);
108    }
109    return s2;
110  }
111
112}}} // of namespaces regression, yat, and theplu
Note: See TracBrowser for help on using the repository browser.