source: trunk/yat/classifier/MatrixLookupWeighted.cc @ 1581

Last change on this file since 1581 was 1581, checked in by Peter, 15 years ago

refs #396 - fixing in KernelLookup?

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date ID
File size: 9.3 KB
Line 
1// $Id$
2
3/*
4  Copyright (C) 2006 Jari Häkkinen, Peter Johansson, Markus Ringnér
5  Copyright (C) 2007 Jari Häkkinen, Peter Johansson
6  Copyright (C) 2008 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 "MatrixLookupWeighted.h"
25#include "MatrixLookup.h"
26#include "yat/utility/DataIterator.h"
27#include "yat/utility/Matrix.h"
28#include "yat/utility/MatrixWeighted.h"
29#include "yat/utility/WeightIterator.h"
30
31#include <algorithm>
32#include <cassert>
33#include <fstream>
34#include <vector>
35
36namespace theplu {
37namespace yat {
38namespace classifier {
39
40  MatrixLookupWeighted::MatrixLookupWeighted(const utility::MatrixWeighted& m,
41                                             const utility::Index& rows,
42                                             const utility::Index& columns)
43    : column_index_(columns), row_index_(rows)
44  {
45    utility::Matrix* data = new utility::Matrix(m.rows(), m.columns());
46    utility::Matrix* weight = new utility::Matrix(m.rows(), m.columns());
47    for (size_t i=0; i<m.rows(); ++i)
48      for (size_t j=0; j<m.columns(); ++j) {
49        (*data)(i,j) = m(i,j).data();
50        (*weight)(i,j) = m(i,j).weight();
51      }
52    // smart pointers are taking ownership
53    data_ = MatrixP(data);
54    weights_ = MatrixP(weight);
55  }
56
57
58  MatrixLookupWeighted::MatrixLookupWeighted(const utility::MatrixWeighted& m,
59                                             bool owner)
60    : column_index_(utility::Index(m.columns())), 
61      row_index_(utility::Index(m.rows()))
62  {
63    // Peter, remember to take care of ownership (but for now leave a leak)
64
65    utility::Matrix* data = new utility::Matrix(m.rows(), m.columns());
66    utility::Matrix* weight = new utility::Matrix(m.rows(), m.columns());
67    for (size_t i=0; i<m.rows(); ++i)
68      for (size_t j=0; j<m.columns(); ++j) {
69        (*data)(i,j) = m(i,j).data();
70        (*weight)(i,j) = m(i,j).weight();
71      } 
72  // smart pointers are taking ownership
73    data_ = MatrixP(data);
74    weights_ = MatrixP(weight);
75  }
76
77
78  MatrixLookupWeighted::MatrixLookupWeighted(const utility::Matrix& data, 
79                                             const utility::Matrix& weights,
80                                             const bool own)
81    : data_(MatrixP(&data, own)), weights_(MatrixP(&weights, own)) 
82  {
83    assert(data.rows()==weights.rows());
84    assert(data.columns()==weights.columns());
85    row_index_ = utility::Index(data.rows());
86    column_index_ = utility::Index(data.columns());
87    assert(validate());
88  }
89
90
91  MatrixLookupWeighted::MatrixLookupWeighted(const MatrixLookup& ml)
92    : column_index_(ml.column_index_), data_(ml.data_),
93      row_index_(ml.row_index_), 
94      weights_(MatrixP(new utility::Matrix(data_->rows(),data_->columns(),1.0)))
95  { 
96    assert(validate());
97  }
98 
99
100  MatrixLookupWeighted::MatrixLookupWeighted(const utility::Matrix& data, 
101                                             const utility::Matrix& weights, 
102                                             const utility::Index& row, 
103                                             const utility::Index& col)
104    : column_index_(col), data_(MatrixP(&data, false)), 
105      row_index_(row), weights_(MatrixP(&weights, false))
106  {
107    assert(validate());
108  }
109 
110
111
112  /*
113  MatrixLookupWeighted::MatrixLookupWeighted(const utility::Matrix& data,
114                                             const utility::Matrix& weights,
115                                             const utility::Index& index,
116                                             const bool row)
117    : data_(MatrixP(new utility::Matrix(data), false)),
118      weights_(MatrixP(new utility::Matrix(weights), false))
119  {
120    assert(data.rows()==weights.rows());
121    assert(data.columns()==weights.columns());
122    if (row){
123      row_index_=index;
124      column_index_ = utility::Index(data.columns());
125    }
126    else{
127      column_index_=index;
128      row_index_ = utility::Index(data.rows());
129    }
130    assert(validate());
131  }
132  */
133
134  /*
135  MatrixLookupWeighted::MatrixLookupWeighted(const MatrixLookup& dv,
136                                             const MatrixLookup& wv)
137    : DataLookup2D(dv), data_(dv.data_), weights_(dv.data_)
138  {
139  }
140  */
141
142
143  MatrixLookupWeighted::MatrixLookupWeighted(const MatrixLookupWeighted& other)
144    : column_index_(other.column_index_), data_(other.data_), 
145      row_index_(other.row_index_), weights_(other.weights_)
146  {
147    assert(validate());
148  }
149
150
151
152  MatrixLookupWeighted::MatrixLookupWeighted(const MatrixLookupWeighted& other,
153                                             const utility::Index& row, 
154                                             const utility::Index& col)
155    : column_index_(utility::Index(other.column_index_, col)), 
156      data_(other.data_), row_index_(utility::Index(other.row_index_, row)), 
157      weights_(other.weights_)
158  {
159    assert(validate());
160  }
161 
162
163
164  MatrixLookupWeighted::MatrixLookupWeighted(const MatrixLookupWeighted& other, 
165                                             const utility::Index& index, 
166                                             bool row)
167    : data_(other.data_), 
168      weights_(other.weights_)
169  {
170    if (row){
171      row_index_ = utility::Index(other.row_index_, index);
172      column_index_= other.column_index_;
173    }
174    else{
175      column_index_ = utility::Index(other.column_index_, index);
176      row_index_= other.row_index_;
177    }
178    assert(validate());
179  }
180 
181
182
183  MatrixLookupWeighted::MatrixLookupWeighted(const size_t rows, 
184                                             const size_t columns, 
185                                             const double value,
186                                             const double weight)
187    : column_index_(utility::Index(std::vector<size_t>(rows, 0))),
188      data_(MatrixP(new utility::Matrix(1,1,value))),
189      row_index_(utility::Index(std::vector<size_t>(columns, 0))), 
190      weights_(MatrixP(new utility::Matrix(1,1,weight)))
191  {
192    assert(validate());
193  }
194
195 
196  MatrixLookupWeighted::MatrixLookupWeighted(std::istream& is, char sep)
197    : data_(MatrixP(new utility::Matrix(is,sep)))
198  {
199    column_index_ = utility::Index(data_->columns());
200    row_index_ = utility::Index(data_->rows());
201    utility::Matrix weights;
202    utility::nan(*data_,weights);
203    // Peter, should be possible to avoid this copying
204    weights_= MatrixP(new utility::Matrix(weights));
205    assert(validate());
206  }
207 
208
209  MatrixLookupWeighted::~MatrixLookupWeighted(void)
210  {
211  }
212
213
214
215  MatrixLookupWeighted::const_iterator MatrixLookupWeighted::begin(void) const
216  {
217    return const_iterator(const_iterator::iterator_type(*this, 0, 0), 1);
218  }
219
220
221  MatrixLookupWeighted::const_column_iterator
222  MatrixLookupWeighted::begin_column(size_t i) const
223  {
224    return const_column_iterator(const_column_iterator::iterator_type(*this,0, 
225                                                                      i),
226                                 columns());
227  }
228
229
230  MatrixLookupWeighted::const_iterator
231  MatrixLookupWeighted::begin_row(size_t i) const
232  {
233    return const_row_iterator(const_row_iterator::iterator_type(*this,i, 0), 1);
234  }
235
236
237  size_t MatrixLookupWeighted::columns(void) const
238  {
239    return column_index_.size();
240  }
241
242
243  double MatrixLookupWeighted::data(size_t row, size_t column) const
244  {
245    assert(row<rows());
246    assert(column<columns());
247    assert(row_index_[row]<data_->rows());
248    assert(column_index_[column]<data_->columns());
249    return (*data_)(row_index_[row], column_index_[column]);
250  }
251
252
253
254  MatrixLookupWeighted::const_iterator MatrixLookupWeighted::end(void) const
255  {
256    return const_iterator(const_iterator::iterator_type(*this, rows(), 0), 1);
257  }
258
259
260  MatrixLookupWeighted::const_iterator
261  MatrixLookupWeighted::end_column(size_t i) const
262  {
263    return const_column_iterator(const_column_iterator::iterator_type(*this, 
264                                                                      rows(),i),
265                                 columns());
266  }
267
268
269  MatrixLookupWeighted::const_row_iterator
270  MatrixLookupWeighted::end_row(size_t i) const
271  {
272    return const_row_iterator(const_row_iterator::iterator_type(*this,i+1,0),1);
273  }
274
275
276  size_t MatrixLookupWeighted::rows(void) const
277  {
278    return row_index_.size();
279  }
280
281
282  bool MatrixLookupWeighted::validate(void) const
283  {
284    for (size_t i=0; i<row_index_.size(); ++i)
285      if (row_index_[i]>=data_->rows()) {
286        std::cerr << i << " " << row_index_[i] << std::endl;
287        return false;
288      }
289    for (size_t i=0; i<column_index_.size(); ++i)
290      if (column_index_[i]>=data_->columns())
291        return false;
292    for (size_t i=0; i<row_index_.size(); ++i)
293      if (row_index_[i]>=weights_->rows())
294        return false;
295    for (size_t i=0; i<column_index_.size(); ++i)
296      if (column_index_[i]>=weights_->columns())
297        return false;
298    return true;
299  }
300
301
302  double MatrixLookupWeighted::weight(size_t row, size_t column) const
303  {
304    assert(row<rows());
305    assert(column<columns());
306    assert(row_index_[row]<weights_->rows());
307    assert(column_index_[column]<weights_->columns());
308    return (*weights_)(row_index_[row], column_index_[column]);
309  }
310
311
312
313  bool MatrixLookupWeighted::weighted(void) const 
314  {
315    return true;
316  }
317
318
319
320  MatrixLookupWeighted::const_reference
321  MatrixLookupWeighted::operator()(const size_t row, const size_t column) const
322  { 
323    return utility::DataWeight(data(row, column), weight(row,column));
324  }
325
326
327
328  const MatrixLookupWeighted& MatrixLookupWeighted::operator=
329  (const MatrixLookupWeighted& other)
330  {
331    if (this!=&other){
332      column_index_=other.column_index_;
333      row_index_=other.row_index_;
334      data_ = other.data_;
335      weights_ = other.weights_;
336    }
337    assert(validate());
338    return *this;
339  }
340
341
342}}} // of namespace classifier, yat, and theplu
Note: See TracBrowser for help on using the repository browser.