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

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

Addresses #436.

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 2.9 KB
RevLine 
[386]1// $Id: MultiDimensional.cc 1486 2008-09-09 21:17:19Z jari $
2
[675]3/*
[831]4  Copyright (C) 2005 Jari Häkkinen
[1275]5  Copyright (C) 2006, 2007, 2008 Jari Häkkinen, Peter Johansson
[386]6
[1437]7  This file is part of the yat library, http://dev.thep.lu.se/yat
[675]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
[1486]11  published by the Free Software Foundation; either version 3 of the
[675]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 this program; if not, write to the Free Software
21  Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA
22  02111-1307, USA.
23*/
24
[680]25#include "MultiDimensional.h"
[750]26#include "yat/utility/Exception.h"
[1121]27#include "yat/utility/Matrix.h"
[1021]28#include "yat/utility/VectorBase.h"
[1120]29#include "yat/utility/Vector.h"
[675]30
[731]31#include <cassert>
32
[386]33namespace theplu {
[680]34namespace yat {
[386]35namespace regression {
36
37
[703]38  MultiDimensional::MultiDimensional(void)
39    : chisquare_(0), work_(NULL)
40  {
41  }
42
43
44  MultiDimensional::~MultiDimensional(void)
45  {
46    if (work_)
47      gsl_multifit_linear_free(work_);
48  }
49
50
[1121]51  const utility::Matrix& MultiDimensional::covariance(void) const
[726]52  {
53    return covariance_;
54  }
55
56
[1121]57  void MultiDimensional::fit(const utility::Matrix& x, 
[1021]58                             const utility::VectorBase& y)
[386]59  {
[731]60    assert(x.rows()==y.size());
[1098]61    covariance_.resize(x.columns(),x.columns());
[1120]62    fit_parameters_ = utility::Vector(x.columns());
[386]63    if (work_)
64      gsl_multifit_linear_free(work_);
[750]65    if (!(work_=gsl_multifit_linear_alloc(x.rows(),fit_parameters_.size())))
66      throw utility::GSL_error("MultiDimensional::fit failed to allocate memory");
67
68    int status = gsl_multifit_linear(x.gsl_matrix_p(), y.gsl_vector_p(),
69                                     fit_parameters_.gsl_vector_p(),
70                                     covariance_.gsl_matrix_p(), &chisquare_,
71                                     work_);
72    if (status)
73      throw utility::GSL_error(std::string("MultiDimensional::fit",status));
74
[739]75    s2_ = chisquare_/(x.rows()-x.columns());
[386]76  }
77
[1120]78  const utility::Vector& MultiDimensional::fit_parameters(void) const
[713]79  { 
80    return fit_parameters_; 
81  }
[386]82
[713]83
84  double MultiDimensional::chisq(void) const
85  {
86    return chisquare_;
87  }
88
89
[1021]90  double MultiDimensional::predict(const utility::VectorBase& x) const 
[718]91  {
[731]92    assert(x.size()==fit_parameters_.size());
[718]93    return fit_parameters_ * x;
94  }
95
96
[1021]97  double MultiDimensional::prediction_error2(const utility::VectorBase& x) const
[586]98  {
[739]99    return standard_error2(x) + s2_;
[586]100  }
101
102
[1021]103  double MultiDimensional::standard_error2(const utility::VectorBase& x) const
[586]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    }
[727]111    return s2;
[586]112  }
113
[681]114}}} // of namespaces regression, yat, and theplu
Note: See TracBrowser for help on using the repository browser.