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

Last change on this file since 2564 was 2564, checked in by Peter, 10 years ago

update copyright years

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