yat  0.18.2pre
Matrix.h
1 #ifndef _theplu_yat_utility_matrix_
2 #define _theplu_yat_utility_matrix_
3 
4 // $Id: Matrix.h 3999 2020-10-08 23:22:32Z peter $
5 
6 /*
7  Copyright (C) 2003 Daniel Dalevi, Peter Johansson
8  Copyright (C) 2004 Jari Häkkinen, Peter Johansson
9  Copyright (C) 2005, 2006 Jari Häkkinen, Peter Johansson, Markus Ringnér
10  Copyright (C) 2007, 2008, 2009 Jari Häkkinen, Peter Johansson
11  Copyright (C) 2012, 2017, 2018, 2019, 2020 Peter Johansson
12 
13  This file is part of the yat library, http://dev.thep.lu.se/yat
14 
15  The yat library is free software; you can redistribute it and/or
16  modify it under the terms of the GNU General Public License as
17  published by the Free Software Foundation; either version 3 of the
18  License, or (at your option) any later version.
19 
20  The yat library is distributed in the hope that it will be useful,
21  but WITHOUT ANY WARRANTY; without even the implied warranty of
22  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
23  General Public License for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with yat. If not, see <http://www.gnu.org/licenses/>.
27 */
28 
29 #include <cassert>
30 
31 #include "config_public.h"
32 
33 // include these operations as they are expected when using Matrix
34 #include "BLAS_level2.h"
35 #include "BLAS_level3.h"
36 
37 #include "Exception.h"
38 #include "MatrixExpression.h"
39 #include "StrideIterator.h"
40 #include "Vector.h"
41 #include "VectorConstView.h"
42 #include "VectorExpression.h"
43 #include "VectorView.h"
44 #include "yat_assert.h"
45 
46 #include <gsl/gsl_matrix.h>
47 
48 #include <cstddef> // size_t
49 #include <iosfwd>
50 
51 namespace theplu {
52 namespace yat {
53 namespace utility {
54 
55  class VectorBase;
56 
74  class Matrix : public BasicMatrix<Matrix>
75  {
76  public:
82  typedef double value_type;
83 
89  typedef double& reference;
90 
96  typedef const double& const_reference;
97 
102 
107 
112 
117 
122 
127 
134  Matrix(void);
135 
144  Matrix(const size_t& r, const size_t& c, double init_value=0);
145 
152  Matrix(const Matrix&);
153 
172  template<class T>
173  Matrix(const MatrixExpression<T>& other);
174 
180  Matrix(Matrix&&) noexcept;
181 
197  template<class T>
199  : blas_result_(NULL), m_(NULL)
200  {
201  other.move(m_);
202  }
203 
204 
235  explicit Matrix(std::istream &, char sep='\0');
236 
240  ~Matrix(void);
241 
245  void all(const double value);
246 
253  iterator begin(void);
254 
261  const_iterator begin(void) const;
262 
268  iterator begin_column(size_t i);
269 
275  const_iterator begin_column(size_t i) const;
276 
282  iterator begin_row(size_t i);
283 
289  const_iterator begin_row(size_t i) const;
290 
294  VectorView column_view(size_t i);
295 
299  const VectorConstView column_const_view(size_t) const;
300 
304  size_t columns(void) const;
305 
313  void div(const Matrix& b);
314 
318  iterator end(void);
319 
323  const_iterator end(void) const;
324 
328  iterator end_column(size_t i);
329 
333  const_iterator end_column(size_t i) const;
334 
338  iterator end_row(size_t i);
339 
343  const_iterator end_row(size_t i) const;
344 
354  bool equal(const Matrix&, const double precision=0) const;
355 
359  const gsl_matrix* gsl_matrix_p(void) const;
360 
364  gsl_matrix* gsl_matrix_p(void);
365 
373  void mul(const Matrix& b);
374 
385  void resize(size_t r, size_t c, double init_value=0);
386 
390  size_t rows(void) const;
391 
395  VectorView row_view(size_t);
396 
400  const VectorConstView row_const_view(size_t) const;
401 
407  void swap_columns(const size_t i,const size_t j);
408 
417  void swap_rowcol(const size_t i,const size_t j);
418 
424  void swap_rows(const size_t i, const size_t j);
425 
435  void transpose(void);
436 
446  double& operator()(size_t row,size_t column);
447 
458  const double& operator()(size_t row,size_t column) const;
459 
471  bool operator==(const Matrix& other) const;
472 
484  bool operator!=(const Matrix& other) const;
485 
494  const Matrix& operator=(const Matrix& other);
495 
514  template<class T>
515  Matrix& operator=(const MatrixExpression<T>& rhs);
516 
517 
523  Matrix& operator=(Matrix&& other);
524 
543  template<class T>
545  {
546  // This function steals data from rhs into this, and give the old
547  // m_ in return to rhs so destructor of rhs takes care of
548  // deallocation.
549 
550  rhs.move(m_);
551  return *this;
552  }
553 
554 
566  const Matrix& operator+=(const Matrix& b);
567 
573  template<class T>
575  {
576  *this = *this + rhs;
577  return *this;
578  }
579 
586  const Matrix& operator+=(const double d);
587 
599  const Matrix& operator-=(const Matrix&);
600 
606  template<class T>
608  {
609  *this = *this - rhs;
610  return *this;
611  }
612 
619  const Matrix& operator-=(const double d);
620 
630  const Matrix& operator*=(const Matrix&);
631 
640  const Matrix& operator*=(double d);
641 
642  private:
643 
655  gsl_matrix* create_gsl_matrix_copy(void) const;
656 
662  void delete_allocated_memory(void);
663 
664  // blas_result_ is used to temporarily store result in BLAS calls
665  // and when not NULL it should always have the same size as
666  // m_. Memory is not allocated for blas_result_ until it is
667  // needed.
668  gsl_matrix* blas_result_;
669  gsl_matrix* m_;
670  };
671 
679  void inverse_svd(const Matrix& input, Matrix& result);
680 
689  bool isnull(const Matrix&);
690 
698  double max(const Matrix&);
699 
707  double min(const Matrix&);
708 
719  void minmax_index(const Matrix&,
720  std::pair<size_t,size_t>& min,
721  std::pair<size_t,size_t>& max);
722 
738  bool nan(const Matrix& templat, Matrix& flag);
739 
753  void swap(Matrix&, Matrix&);
754 
760  std::ostream& operator<< (std::ostream& s, const Matrix&);
761 
771  double trace(const Matrix&);
772 
774  template<class T>
776  : blas_result_(NULL), m_(NULL)
777  {
778  // In principle, we could implement this as a move constructor
779  // and steal the data (rather than copy), but it's a bit hackish
780  // to workaround the constness, so users who are eager for speed
781  // should enable rvalues and get access to the faster move
782  // constructor below.
783 
784  try {
785  detail::copy(m_, other.gsl_matrix_p());
786  }
787  catch (GSL_error& e) {
788  detail::deallocate(m_);
789  throw e;
790  }
791  }
792 
793 
794  template<class T>
796  {
797  // Be careful because *this might be embedded in rhs.
798  rhs.get(blas_result_);
799  std::swap(m_, blas_result_);
800  return *this;
801  }
802 
803 
804 }}} // of namespace utility, yat, and theplu
805 
806 #endif
void mul(const Matrix &b)
const double & const_reference
Definition: Matrix.h:96
iterator end_row(size_t i)
Class for errors reported from underlying GSL calls.
Definition: Exception.h:87
The Department of Theoretical Physics namespace as we define it.
bool operator!=(const Matrix &other) const
Comparison operator. Takes squared time.
StrideIterator< const double * > const_iterator
Definition: Matrix.h:106
This is the yat interface to gsl_vector_view.
Definition: VectorView.h:79
bool operator==(const Matrix &other) const
Comparison operator. Takes squared time.
const gsl_matrix * gsl_matrix_p(void) const
Definition: MatrixExpression.h:138
void resize(size_t r, size_t c, double init_value=0)
Resize Matrix.
iterator begin_row(size_t i)
StrideIterator< const double * > const_row_iterator
Definition: Matrix.h:126
const VectorConstView column_const_view(size_t) const
double & operator()(size_t row, size_t column)
Element access operator.
bool equal(const Matrix &, const double precision=0) const
Check whether matrices are equal within a user defined precision, set by precision.
StrideIterator< double * > row_iterator
Definition: Matrix.h:121
T max(const T &a, const T &b, const T &c)
Definition: stl_utility.h:699
double & reference
Definition: Matrix.h:89
Matrix(void)
The default constructor.
void all(const double value)
size_t rows(void) const
Matrix & operator-=(const MatrixExpression< T > &rhs)
Subtraction and assign operator.
Definition: Matrix.h:607
Matrix & operator+=(const MatrixExpression< T > &rhs)
Addition and assign operator.
Definition: Matrix.h:574
Read-only view.
Definition: VectorConstView.h:56
iterator begin_column(size_t i)
double value_type
Definition: Matrix.h:82
Matrix & operator=(MatrixExpression< T > &&rhs)
Definition: Matrix.h:544
StrideIterator< const double * > const_column_iterator
Definition: Matrix.h:116
~Matrix(void)
The destructor.
const Matrix & operator-=(const Matrix &)
Subtract and assign operator.
const VectorConstView row_const_view(size_t) const
VectorView column_view(size_t i)
An expression that can be converted to a Matrix.
Definition: MatrixExpression.h:46
void div(const Matrix &b)
VectorView row_view(size_t)
Interface to GSL matrix.
Definition: Matrix.h:74
const gsl_matrix * gsl_matrix_p(void) const
void get(gsl_matrix *&m) const
Definition: MatrixExpression.h:110
Definition: BasicMatrix.h:38
const Matrix & operator=(const Matrix &other)
The assignment operator.
StrideIterator< double * > iterator
Definition: Matrix.h:101
StrideIterator< double * > column_iterator
Definition: Matrix.h:111
const Matrix & operator+=(const Matrix &b)
Add and assign operator.
void swap_rows(const size_t i, const size_t j)
Swap rows i and j.
size_t columns(void) const
const Matrix & operator*=(const Matrix &)
Multiply and assignment operator.
iterator end_column(size_t i)
void swap_columns(const size_t i, const size_t j)
Swap columns i and j.
void transpose(void)
Transpose the matrix.
Adaptor using a stride on underlying iterator.
Definition: StrideIterator.h:50
void swap_rowcol(const size_t i, const size_t j)
Swap row i and column j.

Generated on Tue Sep 7 2021 17:32:33 for yat by  doxygen 1.8.14