yat  0.21pre
MatrixBase.h
1 #ifndef _theplu_yat_utility_matrix_base_
2 #define _theplu_yat_utility_matrix_base_
3 
4 // $Id: MatrixBase.h 4252 2022-11-18 02:54:04Z 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, 2021, 2022 Peter Johansson
12 
13  This file is part of the yat library, https://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 <https://www.gnu.org/licenses/>.
27 */
28 
29 #include "config_public.h"
30 
31 // include these operations as they are expected when using Matrix
32 #include "BLAS_level2.h"
33 #include "BLAS_level3.h"
34 
35 #include "Container2DIterator.h"
36 #include "MatrixExpression.h"
37 #include "Vector.h"
38 #include "VectorConstView.h"
39 #include "VectorExpression.h"
40 #include "VectorView.h"
41 #include "yat_assert.h"
42 
43 #include <gsl/gsl_matrix.h>
44 
45 #include <cstddef> // size_t
46 
47 namespace theplu {
48 namespace yat {
49 namespace utility {
50 
54  class MatrixBase : public BasicMatrix<MatrixBase>
55  {
56  public:
62  typedef double value_type;
63 
69  typedef double& reference;
70 
76  typedef const double& const_reference;
77 
83 
88 
93 
100  MatrixBase(void) = default;
101 
105  virtual ~MatrixBase(void);
106 
113  const_iterator begin(void) const;
114 
120  const_column_iterator begin_column(size_t i) const;
121 
127  const_row_iterator begin_row(size_t i) const;
128 
132  const VectorConstView column_const_view(size_t) const;
133 
137  size_t columns(void) const;
138 
142  const_iterator end(void) const;
143 
147  const_column_iterator end_column(size_t i) const;
148 
152  const_row_iterator end_row(size_t i) const;
153 
163  bool equal(const MatrixBase&, const double precision=0) const;
164 
168  virtual const gsl_matrix* gsl_matrix_p(void) const=0;
169 
173  size_t rows(void) const;
174 
178  const VectorConstView row_const_view(size_t) const;
179 
190  const double& operator()(size_t row,size_t column) const;
191 
203  bool operator==(const MatrixBase& other) const;
204 
216  bool operator!=(const MatrixBase& other) const;
217 
218  protected:
230  gsl_matrix* create_gsl_matrix_copy(void) const;
231 
232  };
233 
242  bool isnull(const MatrixBase&);
243 
251  double max(const MatrixBase&);
252 
260  double min(const MatrixBase&);
261 
274  void minmax_index(const MatrixBase& M,
275  std::pair<size_t,size_t>& min,
276  std::pair<size_t,size_t>& max);
277 
283  std::ostream& operator<< (std::ostream& s, const MatrixBase&);
284 
294  double trace(const MatrixBase&);
295 
296 }}} // of namespace utility, yat, and theplu
297 
298 #endif
const_iterator begin(void) const
Definition: MatrixBase.h:54
const_row_iterator begin_row(size_t i) const
Container2DIterator< const MatrixBase, const double, const double & > const_iterator
Definition: MatrixBase.h:82
The Department of Theoretical Physics namespace as we define it.
const_row_iterator end_row(size_t i) const
MatrixBase(void)=default
The default constructor.
gsl_matrix * create_gsl_matrix_copy(void) const
Create a new copy of the internal GSL matrix.
double & reference
Definition: MatrixBase.h:69
StrideIterator< const double * > const_row_iterator
Definition: MatrixBase.h:92
const_column_iterator begin_column(size_t i) const
bool operator==(const MatrixBase &other) const
Comparison operator. Takes squared time.
T max(const T &a, const T &b, const T &c)
Definition: stl_utility.h:699
Read-only view.
Definition: VectorConstView.h:56
Iterator for a Container2D.
Definition: Container2DIterator.h:62
const_iterator end(void) const
double value_type
Definition: MatrixBase.h:62
bool operator!=(const MatrixBase &other) const
Comparison operator. Takes squared time.
virtual ~MatrixBase(void)
The destructor.
bool equal(const MatrixBase &, const double precision=0) const
Check whether matrices are equal within a user defined precision, set by precision.
const VectorConstView column_const_view(size_t) const
Definition: BasicMatrix.h:38
virtual const gsl_matrix * gsl_matrix_p(void) const =0
const double & operator()(size_t row, size_t column) const
Element access operator.
const VectorConstView row_const_view(size_t) const
StrideIterator< const double * > const_column_iterator
Definition: MatrixBase.h:87
const double & const_reference
Definition: MatrixBase.h:76
const_column_iterator end_column(size_t i) const
Adaptor using a stride on underlying iterator.
Definition: StrideIterator.h:50

Generated on Wed Jan 25 2023 03:34:29 for yat by  doxygen 1.8.14