1 | // $Id: Matrix.cc 1681 2008-12-29 20:47:47Z peter $ |
---|
2 | |
---|
3 | /* |
---|
4 | Copyright (C) 2003 Daniel Dalevi, Peter Johansson |
---|
5 | Copyright (C) 2004 Jari Häkkinen, Peter Johansson |
---|
6 | Copyright (C) 2005, 2006 Jari Häkkinen, Peter Johansson, Markus Ringnér |
---|
7 | Copyright (C) 2007, 2008 Jari Häkkinen, Peter Johansson |
---|
8 | |
---|
9 | This file is part of the yat library, http://dev.thep.lu.se/yat |
---|
10 | |
---|
11 | The yat library is free software; you can redistribute it and/or |
---|
12 | modify it under the terms of the GNU General Public License as |
---|
13 | published by the Free Software Foundation; either version 3 of the |
---|
14 | License, or (at your option) any later version. |
---|
15 | |
---|
16 | The yat library is distributed in the hope that it will be useful, |
---|
17 | but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
18 | MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
---|
19 | General Public License for more details. |
---|
20 | |
---|
21 | You should have received a copy of the GNU General Public License |
---|
22 | along with yat. If not, see <http://www.gnu.org/licenses/>. |
---|
23 | */ |
---|
24 | |
---|
25 | #include "Matrix.h" |
---|
26 | #include "Vector.h" |
---|
27 | #include "VectorBase.h" |
---|
28 | #include "VectorConstView.h" |
---|
29 | #include "VectorView.h" |
---|
30 | #include "utility.h" |
---|
31 | |
---|
32 | #include <cassert> |
---|
33 | #include <cmath> |
---|
34 | #include <sstream> |
---|
35 | #include <vector> |
---|
36 | |
---|
37 | #include <gsl/gsl_blas.h> |
---|
38 | |
---|
39 | namespace theplu { |
---|
40 | namespace yat { |
---|
41 | namespace utility { |
---|
42 | |
---|
43 | |
---|
44 | Matrix::Matrix(void) |
---|
45 | : blas_result_(NULL), m_(NULL) |
---|
46 | { |
---|
47 | } |
---|
48 | |
---|
49 | |
---|
50 | Matrix::Matrix(const size_t& r, const size_t& c, double init_value) |
---|
51 | : blas_result_(NULL), m_(NULL) |
---|
52 | { |
---|
53 | resize(r,c,init_value); |
---|
54 | } |
---|
55 | |
---|
56 | |
---|
57 | Matrix::Matrix(const Matrix& o) |
---|
58 | : blas_result_(NULL), m_(o.create_gsl_matrix_copy()) |
---|
59 | { |
---|
60 | } |
---|
61 | |
---|
62 | |
---|
63 | // Constructor that gets data from istream |
---|
64 | Matrix::Matrix(std::istream& is, char sep) |
---|
65 | throw (utility::IO_error,std::exception) |
---|
66 | : blas_result_(NULL) |
---|
67 | { |
---|
68 | if (!is.good()) |
---|
69 | throw utility::IO_error("Matrix: istream is not good"); |
---|
70 | |
---|
71 | // read the data file and store in stl vectors (dynamically |
---|
72 | // expandable) |
---|
73 | std::vector<std::vector<double> > data_matrix; |
---|
74 | unsigned int nof_columns=0; |
---|
75 | unsigned int nof_rows = 0; |
---|
76 | std::string line; |
---|
77 | while(getline(is, line, '\n')){ |
---|
78 | // Ignoring empty lines |
---|
79 | if (!line.size()) { |
---|
80 | continue; |
---|
81 | } |
---|
82 | nof_rows++; |
---|
83 | data_matrix.resize(data_matrix.size()+1); |
---|
84 | std::vector<double>& v=data_matrix.back(); |
---|
85 | std::string element; |
---|
86 | std::stringstream ss(line); |
---|
87 | |
---|
88 | bool ok=true; |
---|
89 | while(ok) { |
---|
90 | if(sep=='\0') |
---|
91 | ok=(ss>>element); |
---|
92 | else |
---|
93 | ok=getline(ss, element, sep); |
---|
94 | if(!ok) |
---|
95 | break; |
---|
96 | |
---|
97 | if (!element.size()) |
---|
98 | v.push_back(std::numeric_limits<double>::quiet_NaN()); |
---|
99 | else { |
---|
100 | try { |
---|
101 | v.push_back(convert<double>(element)); |
---|
102 | } |
---|
103 | catch (std::runtime_error& e) { |
---|
104 | std::stringstream ss(e.what()); |
---|
105 | ss << "\nMatrix.cc: " << element |
---|
106 | << " is not accepted as a matrix element\n"; |
---|
107 | throw IO_error(ss.str()); |
---|
108 | } |
---|
109 | } |
---|
110 | } |
---|
111 | if(sep!='\0' && line[line.size()-1]==sep) // add NaN for final separator |
---|
112 | v.push_back(std::numeric_limits<double>::quiet_NaN()); |
---|
113 | if (!nof_columns) |
---|
114 | nof_columns=v.size(); |
---|
115 | else if (v.size()!=nof_columns) { |
---|
116 | std::ostringstream s; |
---|
117 | s << "Matrix::Matrix(std::istream&, char) data file error: " |
---|
118 | << "line " << nof_rows << " has " << v.size() |
---|
119 | << " columns; expected " << nof_columns << " columns."; |
---|
120 | throw utility::IO_error(s.str()); |
---|
121 | } |
---|
122 | } |
---|
123 | |
---|
124 | // manipulate the state of the stream to be good |
---|
125 | is.clear(std::ios::goodbit); |
---|
126 | |
---|
127 | // if stream was empty, create nothing |
---|
128 | if (!nof_columns || !nof_rows) |
---|
129 | return; |
---|
130 | |
---|
131 | // convert the data to a gsl matrix |
---|
132 | m_ = gsl_matrix_alloc ( nof_rows, nof_columns ); |
---|
133 | if (!m_) |
---|
134 | throw utility::GSL_error("Matrix::Matrix failed to allocate memory"); |
---|
135 | |
---|
136 | // if gsl error handler disabled, out of bounds index will not |
---|
137 | // abort the program. |
---|
138 | for(size_t i=0;i<nof_rows;i++) |
---|
139 | for(size_t j=0;j<nof_columns;j++) |
---|
140 | gsl_matrix_set( m_, i, j, data_matrix[i][j] ); |
---|
141 | } |
---|
142 | |
---|
143 | |
---|
144 | Matrix::~Matrix(void) |
---|
145 | { |
---|
146 | delete_allocated_memory(); |
---|
147 | if (blas_result_) |
---|
148 | gsl_matrix_free(blas_result_); |
---|
149 | blas_result_=NULL; |
---|
150 | } |
---|
151 | |
---|
152 | |
---|
153 | void Matrix::all(const double value) |
---|
154 | { |
---|
155 | assert(m_); |
---|
156 | gsl_matrix_set_all(m_, value); |
---|
157 | } |
---|
158 | |
---|
159 | |
---|
160 | Matrix::iterator Matrix::begin(void) |
---|
161 | { |
---|
162 | return iterator(&(*this)(0,0), 1); |
---|
163 | } |
---|
164 | |
---|
165 | |
---|
166 | Matrix::const_iterator Matrix::begin(void) const |
---|
167 | { |
---|
168 | return const_iterator(&(*this)(0,0), 1); |
---|
169 | } |
---|
170 | |
---|
171 | |
---|
172 | Matrix::column_iterator Matrix::begin_column(size_t i) |
---|
173 | { |
---|
174 | return iterator(&(*this)(0,i), this->columns()); |
---|
175 | } |
---|
176 | |
---|
177 | |
---|
178 | Matrix::const_column_iterator Matrix::begin_column(size_t i) const |
---|
179 | { |
---|
180 | return const_iterator(&(*this)(0,i), this->columns()); |
---|
181 | } |
---|
182 | |
---|
183 | |
---|
184 | Matrix::row_iterator Matrix::begin_row(size_t i) |
---|
185 | { |
---|
186 | return iterator(&(*this)(i,0), 1); |
---|
187 | } |
---|
188 | |
---|
189 | |
---|
190 | Matrix::const_row_iterator Matrix::begin_row(size_t i) const |
---|
191 | { |
---|
192 | return const_iterator(&(*this)(i,0), 1); |
---|
193 | } |
---|
194 | |
---|
195 | |
---|
196 | VectorView Matrix::column_view(size_t col) |
---|
197 | { |
---|
198 | VectorView res(*this, col, false); |
---|
199 | return res; |
---|
200 | } |
---|
201 | |
---|
202 | |
---|
203 | const VectorConstView Matrix::column_const_view(size_t col) const |
---|
204 | { |
---|
205 | return VectorConstView(*this, col, false); |
---|
206 | } |
---|
207 | |
---|
208 | |
---|
209 | size_t Matrix::columns(void) const |
---|
210 | { |
---|
211 | return (m_ ? m_->size2 : 0); |
---|
212 | } |
---|
213 | |
---|
214 | |
---|
215 | gsl_matrix* Matrix::create_gsl_matrix_copy(void) const |
---|
216 | { |
---|
217 | if (!m_) |
---|
218 | return NULL; |
---|
219 | gsl_matrix* m = gsl_matrix_alloc(rows(),columns()); |
---|
220 | if (!m) |
---|
221 | throw utility::GSL_error("Matrix::create_gsl_matrix_copy failed to allocate memory"); |
---|
222 | if (gsl_matrix_memcpy(m,m_)) |
---|
223 | throw utility::GSL_error("Matrix::create_gsl_matrix_copy dimension mis-match"); |
---|
224 | return m; |
---|
225 | } |
---|
226 | |
---|
227 | |
---|
228 | void Matrix::delete_allocated_memory(void) |
---|
229 | { |
---|
230 | if (m_) |
---|
231 | gsl_matrix_free(m_); |
---|
232 | m_=NULL; |
---|
233 | } |
---|
234 | |
---|
235 | |
---|
236 | void Matrix::div(const Matrix& other) |
---|
237 | { |
---|
238 | assert(m_); |
---|
239 | int status=gsl_matrix_div_elements(m_, other.gsl_matrix_p()); |
---|
240 | if (status) |
---|
241 | throw utility::GSL_error(std::string("matrix::div_elements",status)); |
---|
242 | } |
---|
243 | |
---|
244 | |
---|
245 | Matrix::iterator Matrix::end(void) |
---|
246 | { |
---|
247 | return iterator(&(*this)(0,0)+rows()*columns(), 1); |
---|
248 | } |
---|
249 | |
---|
250 | |
---|
251 | Matrix::const_iterator Matrix::end(void) const |
---|
252 | { |
---|
253 | return const_iterator(&(*this)(0,0)+rows()*columns(), 1); |
---|
254 | } |
---|
255 | |
---|
256 | |
---|
257 | Matrix::column_iterator Matrix::end_column(size_t i) |
---|
258 | { |
---|
259 | return column_iterator(&(*this)(0,i)+rows()*columns(), this->columns()); |
---|
260 | } |
---|
261 | |
---|
262 | |
---|
263 | Matrix::const_column_iterator Matrix::end_column(size_t i) const |
---|
264 | { |
---|
265 | return const_column_iterator(&(*this)(0,i)+rows()*columns(),this->columns()); |
---|
266 | } |
---|
267 | |
---|
268 | |
---|
269 | Matrix::row_iterator Matrix::end_row(size_t i) |
---|
270 | { |
---|
271 | return row_iterator(&(*this)(i,0)+columns(), 1); |
---|
272 | } |
---|
273 | |
---|
274 | |
---|
275 | Matrix::const_row_iterator Matrix::end_row(size_t i) const |
---|
276 | { |
---|
277 | return const_row_iterator(&(*this)(i,0)+columns(), 1); |
---|
278 | } |
---|
279 | |
---|
280 | |
---|
281 | bool Matrix::equal(const Matrix& other, const double d) const |
---|
282 | { |
---|
283 | if (this==&other) |
---|
284 | return true; |
---|
285 | if (columns()!=other.columns() || rows()!=other.rows()) |
---|
286 | return false; |
---|
287 | for (size_t i=0; i<rows(); i++) |
---|
288 | for (size_t j=0; j<columns(); j++) |
---|
289 | // The two last condition checks are needed for NaN detection |
---|
290 | if (fabs( (*this)(i,j)-other(i,j) ) > d || |
---|
291 | (*this)(i,j)!=(*this)(i,j) || other(i,j)!=other(i,j)) |
---|
292 | return false; |
---|
293 | return true; |
---|
294 | } |
---|
295 | |
---|
296 | |
---|
297 | const gsl_matrix* Matrix::gsl_matrix_p(void) const |
---|
298 | { |
---|
299 | return m_; |
---|
300 | } |
---|
301 | |
---|
302 | |
---|
303 | gsl_matrix* Matrix::gsl_matrix_p(void) |
---|
304 | { |
---|
305 | return m_; |
---|
306 | } |
---|
307 | |
---|
308 | |
---|
309 | void Matrix::mul(const Matrix& other) |
---|
310 | { |
---|
311 | assert(m_); |
---|
312 | int status=gsl_matrix_mul_elements(m_, other.gsl_matrix_p()); |
---|
313 | if (status) |
---|
314 | throw utility::GSL_error(std::string("Matrix::mul_elements",status)); |
---|
315 | } |
---|
316 | |
---|
317 | |
---|
318 | void Matrix::resize(size_t r, size_t c, double init_value) |
---|
319 | { |
---|
320 | delete_allocated_memory(); |
---|
321 | |
---|
322 | assert( (r&&c) || (!r&&!c) ); |
---|
323 | if (r && c){ |
---|
324 | m_ = gsl_matrix_alloc(r,c); |
---|
325 | if (!m_) |
---|
326 | throw utility::GSL_error("Matrix::Matrix failed to allocate memory"); |
---|
327 | all(init_value); |
---|
328 | } |
---|
329 | |
---|
330 | // no need to delete blas_result_ if the number of rows fit, it |
---|
331 | // may be useful later. |
---|
332 | if (blas_result_ && (blas_result_->size1!=rows())) { |
---|
333 | gsl_matrix_free(blas_result_); |
---|
334 | blas_result_=NULL; |
---|
335 | } |
---|
336 | } |
---|
337 | |
---|
338 | |
---|
339 | size_t Matrix::rows(void) const |
---|
340 | { |
---|
341 | return (m_ ? m_->size1 : 0); |
---|
342 | } |
---|
343 | |
---|
344 | |
---|
345 | const VectorConstView Matrix::row_const_view(size_t col) const |
---|
346 | { |
---|
347 | return VectorConstView(*this, col, true); |
---|
348 | } |
---|
349 | |
---|
350 | |
---|
351 | VectorView Matrix::row_view(size_t row) |
---|
352 | { |
---|
353 | VectorView res(*this, row, true); |
---|
354 | return res; |
---|
355 | } |
---|
356 | |
---|
357 | |
---|
358 | void Matrix::swap_columns(const size_t i, const size_t j) |
---|
359 | { |
---|
360 | assert(m_); |
---|
361 | int status=gsl_matrix_swap_columns(m_, i, j); |
---|
362 | if (status) |
---|
363 | throw utility::GSL_error(std::string("Matrix::swap_columns",status)); |
---|
364 | } |
---|
365 | |
---|
366 | |
---|
367 | void Matrix::swap_rowcol(const size_t i, const size_t j) |
---|
368 | { |
---|
369 | assert(m_); |
---|
370 | int status=gsl_matrix_swap_rowcol(m_, i, j); |
---|
371 | if (status) |
---|
372 | throw utility::GSL_error(std::string("Matrix::swap_rowcol",status)); |
---|
373 | } |
---|
374 | |
---|
375 | |
---|
376 | void Matrix::swap_rows(const size_t i, const size_t j) |
---|
377 | { |
---|
378 | assert(m_); |
---|
379 | int status=gsl_matrix_swap_rows(m_, i, j); |
---|
380 | if (status) |
---|
381 | throw utility::GSL_error(std::string("Matrix::swap_rows",status)); |
---|
382 | } |
---|
383 | |
---|
384 | |
---|
385 | void Matrix::transpose(void) |
---|
386 | { |
---|
387 | assert(m_); |
---|
388 | if (columns()==rows()) |
---|
389 | gsl_matrix_transpose(m_); // this never fails |
---|
390 | else { |
---|
391 | gsl_matrix* transposed = gsl_matrix_alloc(columns(),rows()); |
---|
392 | if (!transposed) |
---|
393 | throw utility::GSL_error("Matrix::transpose failed to allocate memory"); |
---|
394 | // next line never fails if allocation above succeeded. |
---|
395 | gsl_matrix_transpose_memcpy(transposed,m_); |
---|
396 | gsl_matrix_free(m_); |
---|
397 | m_ = transposed; |
---|
398 | if (blas_result_) { |
---|
399 | gsl_matrix_free(blas_result_); |
---|
400 | blas_result_=NULL; |
---|
401 | } |
---|
402 | } |
---|
403 | } |
---|
404 | |
---|
405 | |
---|
406 | double& Matrix::operator()(size_t row, size_t column) |
---|
407 | { |
---|
408 | assert(m_); |
---|
409 | assert(row<rows()); |
---|
410 | assert(column<columns()); |
---|
411 | double* d=gsl_matrix_ptr(m_, row, column); |
---|
412 | if (!d) |
---|
413 | throw utility::GSL_error("Matrix::operator()",GSL_EINVAL); |
---|
414 | return *d; |
---|
415 | } |
---|
416 | |
---|
417 | |
---|
418 | const double& Matrix::operator()(size_t row, size_t column) const |
---|
419 | { |
---|
420 | assert(row<rows()); |
---|
421 | assert(column<columns()); |
---|
422 | const double* d=gsl_matrix_const_ptr(m_, row, column); |
---|
423 | if (!d) |
---|
424 | throw utility::GSL_error("Matrix::operator()",GSL_EINVAL); |
---|
425 | return *d; |
---|
426 | } |
---|
427 | |
---|
428 | |
---|
429 | bool Matrix::operator==(const Matrix& other) const |
---|
430 | { |
---|
431 | return equal(other); |
---|
432 | } |
---|
433 | |
---|
434 | |
---|
435 | bool Matrix::operator!=(const Matrix& other) const |
---|
436 | { |
---|
437 | return !equal(other); |
---|
438 | } |
---|
439 | |
---|
440 | |
---|
441 | const Matrix& Matrix::operator=( const Matrix& other ) |
---|
442 | { |
---|
443 | if (this!=&other) { |
---|
444 | if (rows()!=other.rows() || columns()!=other.columns()) |
---|
445 | resize(other.rows(), other.columns()); |
---|
446 | |
---|
447 | if (m_) |
---|
448 | if (gsl_matrix_memcpy(m_, other.gsl_matrix_p())) { |
---|
449 | std::string s="Matrix::operator= assignment failed"; |
---|
450 | throw utility::GSL_error(s); |
---|
451 | } |
---|
452 | } |
---|
453 | return *this; |
---|
454 | } |
---|
455 | |
---|
456 | |
---|
457 | const Matrix& Matrix::operator+=(const Matrix& other) |
---|
458 | { |
---|
459 | assert(m_); |
---|
460 | int status=gsl_matrix_add(m_, other.m_); |
---|
461 | if (status) |
---|
462 | throw utility::GSL_error(std::string("Matrix::operator+=", status)); |
---|
463 | return *this; |
---|
464 | } |
---|
465 | |
---|
466 | |
---|
467 | const Matrix& Matrix::operator+=(const double d) |
---|
468 | { |
---|
469 | assert(m_); |
---|
470 | gsl_matrix_add_constant(m_, d); |
---|
471 | return *this; |
---|
472 | } |
---|
473 | |
---|
474 | |
---|
475 | const Matrix& Matrix::operator-=(const Matrix& other) |
---|
476 | { |
---|
477 | assert(m_); |
---|
478 | int status=gsl_matrix_sub(m_, other.m_); |
---|
479 | if (status) |
---|
480 | throw utility::GSL_error(std::string("Matrix::operator-=", status)); |
---|
481 | return *this; |
---|
482 | } |
---|
483 | |
---|
484 | |
---|
485 | const Matrix& Matrix::operator-=(const double d) |
---|
486 | { |
---|
487 | assert(m_); |
---|
488 | gsl_matrix_add_constant(m_, -d); |
---|
489 | return *this; |
---|
490 | } |
---|
491 | |
---|
492 | |
---|
493 | const Matrix& Matrix::operator*=(const Matrix& other) |
---|
494 | { |
---|
495 | assert(m_); |
---|
496 | if ( blas_result_ && ((blas_result_->size1!=rows()) || |
---|
497 | (blas_result_->size2!=other.columns())) ) { |
---|
498 | gsl_matrix_free(blas_result_); |
---|
499 | blas_result_=NULL; |
---|
500 | } |
---|
501 | if (!blas_result_) { |
---|
502 | blas_result_ = gsl_matrix_alloc(rows(),other.columns()); |
---|
503 | if (!blas_result_) |
---|
504 | throw utility::GSL_error("Matrix::operator*= failed to allocate memory"); |
---|
505 | } |
---|
506 | gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, m_, other.m_, 0.0, blas_result_); |
---|
507 | gsl_matrix* tmp=m_; |
---|
508 | m_ = blas_result_; |
---|
509 | blas_result_=tmp; |
---|
510 | return *this; |
---|
511 | } |
---|
512 | |
---|
513 | |
---|
514 | const Matrix& Matrix::operator*=(const double d) |
---|
515 | { |
---|
516 | assert(m_); |
---|
517 | gsl_matrix_scale(m_, d); |
---|
518 | return *this; |
---|
519 | } |
---|
520 | |
---|
521 | |
---|
522 | bool isnull(const Matrix& other) |
---|
523 | { |
---|
524 | return gsl_matrix_isnull(other.gsl_matrix_p()); |
---|
525 | } |
---|
526 | |
---|
527 | |
---|
528 | double max(const Matrix& other) |
---|
529 | { |
---|
530 | return gsl_matrix_max(other.gsl_matrix_p()); |
---|
531 | } |
---|
532 | |
---|
533 | |
---|
534 | double min(const Matrix& other) |
---|
535 | { |
---|
536 | return gsl_matrix_min(other.gsl_matrix_p()); |
---|
537 | } |
---|
538 | |
---|
539 | |
---|
540 | void minmax_index(const Matrix& other, |
---|
541 | std::pair<size_t,size_t>& min, std::pair<size_t,size_t>& max) |
---|
542 | { |
---|
543 | gsl_matrix_minmax_index(other.gsl_matrix_p(), &min.first, &min.second, |
---|
544 | &max.first, &max.second); |
---|
545 | } |
---|
546 | |
---|
547 | |
---|
548 | bool nan(const Matrix& templat, Matrix& flag) |
---|
549 | { |
---|
550 | size_t rows=templat.rows(); |
---|
551 | size_t columns=templat.columns(); |
---|
552 | if (rows!=flag.rows() && columns!=flag.columns()) |
---|
553 | flag.resize(rows,columns,1.0); |
---|
554 | else |
---|
555 | flag.all(1.0); |
---|
556 | bool nan=false; |
---|
557 | for (size_t i=0; i<rows; i++) |
---|
558 | for (size_t j=0; j<columns; j++) |
---|
559 | if (std::isnan(templat(i,j))) { |
---|
560 | flag(i,j)=0; |
---|
561 | nan=true; |
---|
562 | } |
---|
563 | return nan; |
---|
564 | } |
---|
565 | |
---|
566 | |
---|
567 | void swap(Matrix& a, Matrix& b) |
---|
568 | { |
---|
569 | assert(a.gsl_matrix_p()); assert(b.gsl_matrix_p()); |
---|
570 | int status=gsl_matrix_swap(a.gsl_matrix_p(), b.gsl_matrix_p()); |
---|
571 | if (status) |
---|
572 | throw utility::GSL_error(std::string("swap(Matrix&,Matrix&)",status)); |
---|
573 | } |
---|
574 | |
---|
575 | |
---|
576 | std::ostream& operator<<(std::ostream& s, const Matrix& m) |
---|
577 | { |
---|
578 | s.setf(std::ios::dec); |
---|
579 | s.precision(12); |
---|
580 | for(size_t i=0, j=0; i<m.rows(); i++) |
---|
581 | for (j=0; j<m.columns(); j++) { |
---|
582 | s << m(i,j); |
---|
583 | if (j<m.columns()-1) |
---|
584 | s << s.fill(); |
---|
585 | else if (i<m.rows()-1) |
---|
586 | s << "\n"; |
---|
587 | } |
---|
588 | return s; |
---|
589 | } |
---|
590 | |
---|
591 | |
---|
592 | Vector operator*(const Matrix& m, const VectorBase& v) |
---|
593 | { |
---|
594 | utility::Vector res(m.rows()); |
---|
595 | for (size_t i=0; i<res.size(); ++i) |
---|
596 | res(i) = VectorConstView(m,i) * v; |
---|
597 | return res; |
---|
598 | } |
---|
599 | |
---|
600 | |
---|
601 | Vector operator*(const VectorBase& v, const Matrix& m) |
---|
602 | { |
---|
603 | utility::Vector res(m.columns()); |
---|
604 | for (size_t i=0; i<res.size(); ++i) |
---|
605 | res(i) = v * VectorConstView(m,i,false); |
---|
606 | return res; |
---|
607 | } |
---|
608 | |
---|
609 | }}} // of namespace utility, yat and thep |
---|