ampsci
High-precision calculations for one- and two-valence atomic systems
Matrix.hpp
1#pragma once
2#include "qip/Vector.hpp" // for std::vector overloads
3#include <array>
4#include <cassert>
5#include <complex>
6#include <gsl/gsl_blas.h>
7#include <gsl/gsl_complex.h>
8#include <gsl/gsl_complex_math.h>
9#include <gsl/gsl_eigen.h>
10#include <gsl/gsl_linalg.h>
11#include <gsl/gsl_math.h>
12#include <iostream>
13#include <type_traits>
14#include <utility>
15#include <vector>
16
17template <typename T>
18struct is_complex : std::false_type {};
19template <typename T>
20struct is_complex<std::complex<T>> : std::true_type {};
21template <typename T>
22constexpr bool is_complex_v = is_complex<T>::value;
23
24//==============================================================================
25//! Defines Matrix, Vector classes, and linear some algebra functions
26namespace LinAlg {
27
28//! Proved a "view" onto an array
29template <typename T>
30class View;
31
32//! Provides a non-owning view onto a Matrix (read/write, but no resize)
33template <typename T>
34class Matrix_view;
35
36//==============================================================================
37//! Matrix class; row-major
38template <typename T = double>
39class Matrix {
40
41protected:
42 std::size_t m_rows;
43 std::size_t m_cols;
44 std::vector<T> m_data{};
45
46public:
47 //! Default initialiser
48 Matrix() : m_rows(0), m_cols(0) {}
49
50 //! Initialise a blank matrix rows*cols, filled with 0
51 Matrix(std::size_t rows, std::size_t cols)
52 : m_rows(rows), m_cols(cols), m_data(rows * cols) {}
53
54 //! Initialise a matrix rows*cols, filled with 'value'
55 Matrix(std::size_t rows, std::size_t cols, const T &value)
56 : m_rows(rows), m_cols(cols), m_data(rows * cols, value) {}
57
58 //! Initialise a blank square matrix dimension*dimension, filled with 0
59 // excplicit, since don't alow flot->int converions
60 explicit Matrix(std::size_t dimension)
61 : m_rows(dimension), m_cols(dimension), m_data(dimension * dimension) {}
62
63 //! Initialise a matrix from initialiser list. {{},{},{}}. Each row must be
64 //! same length
65 Matrix(std::initializer_list<std::initializer_list<T>> ll)
66 : m_rows(ll.size()), m_cols(ll.begin()->size()), m_data{} {
67 // way to avoid copy?
68 m_data.reserve(m_rows * m_cols);
69 for (auto &l : ll) {
70 m_data.insert(m_data.end(), l.begin(), l.end());
71 }
72 }
73
74 //! Initialise a matrix from single initialiser list. {...}.
75 Matrix(std::size_t rows, std::size_t cols, std::initializer_list<T> l)
76 : m_rows(rows), m_cols(cols), m_data{l} {
77 assert(m_data.size() == rows * cols &&
78 "initializer_list must be rows*cols");
79 }
80
81 //! Initialise a matrix from single initialiser list. {...}.
82 Matrix(std::size_t rows, std::size_t cols, std::vector<T> &&v)
83 : m_rows(rows), m_cols(cols), m_data{std::forward<std::vector<T>>(v)} {
84 assert(m_data.size() == rows * cols &&
85 "initializer_list must be rows*cols");
86 }
87 //! Initialise a matrix from single initialiser list. {...}.
88 Matrix(std::size_t rows, std::size_t cols, const std::vector<T> &v)
89 : m_rows(rows), m_cols(cols), m_data{v} {
90 assert(m_data.size() == rows * cols &&
91 "initializer_list must be rows*cols");
92 }
93
94 //============================================================================
95 //! Resizes matrix to new dimension; all values reset to default
96 void resize(std::size_t rows, std::size_t cols) {
97 m_rows = rows;
98 m_cols = cols;
99 m_data.assign(rows * cols, T{});
100 }
101
102 //! Resizes matrix to new dimension; all values reset to 'value'
103 void resize(std::size_t rows, std::size_t cols, const T &value) {
104 m_rows = rows;
105 m_cols = cols;
106 m_data.assign(rows * cols, value);
107 }
108
109 //============================================================================
110
111 //! Return rows [major index size]
112 std::size_t rows() const { return m_rows; }
113 //! Return columns [minor index size]
114 std::size_t cols() const { return m_cols; }
115 //! Return rows*columns [total array size]
116 std::size_t size() const { return m_data.size(); }
117
118 //! Bool - is empty? (same as size==0)
119 bool empty() const { return m_data.empty(); }
120
121 //! Returns pointer to first element. Note: for std::complex<T>, this is a
122 //! pointer to complex<T>, not T
123 T *data() { return m_data.data(); }
124 //! As above, but const
125 const T *data() const { return m_data.data(); }
126
127 //============================================================================
128
129 //! [] index access (with no range checking). [i][j] returns ith row, jth col
130 const T *operator[](std::size_t i) const { return &(m_data[i * m_cols]); }
131 //! As above, but const
132 T *operator[](std::size_t i) { return &(m_data[i * m_cols]); }
133
134 //! () index access (with range checking). (i,j) returns ith row, jth col
135 T &at(std::size_t row_i, std::size_t col_j) {
136 assert(row_i < m_rows && col_j < m_cols);
137 return m_data[row_i * m_cols + col_j];
138 }
139 //! const () index access (with range checking). (i,j) ith row, jth col
140 T at(std::size_t row_i, std::size_t col_j) const {
141 assert(row_i < m_rows && col_j < m_cols);
142 return m_data[row_i * m_cols + col_j];
143 }
144
145 //! const ref () index access (with range checking). (i,j) ith row, jth col
146 const T &atc(std::size_t row_i, std::size_t col_j) const {
147 assert(row_i < m_rows && col_j < m_cols);
148 return m_data[row_i * m_cols + col_j];
149 }
150
151 //! () index access (with range checking). (i,j) returns ith row, jth col
152 T &operator()(std::size_t i, std::size_t j) { return at(i, j); }
153 //! As above, but const
154 T operator()(std::size_t i, std::size_t j) const { return at(i, j); }
155
156 //============================================================================
157
158 //! iterators for underlying std::vector (entire data)
159 auto begin() { return m_data.begin(); }
160 auto cbegin() const { return m_data.cbegin(); }
161 auto end() { return m_data.end(); }
162 auto cend() const { return m_data.cend(); }
163
164 //! Returns raw c pointer to start of a row
165 // [[deprecated]]
166 const T *row(std::size_t row) const {
167 return m_data.data() + long(row * m_cols);
168 }
169
170 //! Returns a mutable 'View' of a row
171 [[nodiscard]] View<T> row_view(std::size_t row) {
172 return View<T>(this->data(), row * m_cols, m_cols, 1ul);
173 }
174 //! Returns an immutable 'View' of a row
175 [[nodiscard]] View<const T> row_view(std::size_t row) const {
176 return View<const T>(this->data(), row * m_cols, m_cols, 1ul);
177 }
178 //! Returns a mutable 'View' of a column
179 [[nodiscard]] View<T> column_view(std::size_t col) {
180 return View<T>(this->data(), col, m_rows, m_rows);
181 }
182 //! Returns an immutable 'View' of a column
183 [[nodiscard]] View<const T> column_view(std::size_t col) const {
184 return View<const T>(this->data(), col, m_rows, m_rows);
185 }
186
187 //! Returns a mutable view of the entire matrix (no ownership, no resize)
188 [[nodiscard]] Matrix_view<T> matrix_view() {
189 return Matrix_view<T>(m_data.data(), m_rows, m_cols);
190 }
191 //! Returns an immutable view of the entire matrix (no ownership, no resize)
192 [[nodiscard]] Matrix_view<const T> matrix_view() const {
193 return Matrix_view<const T>(m_data.data(), m_rows, m_cols);
194 }
195
196 //============================================================================
197 //! Returns gsl_matrix_view (or _float_view, _complex_view,
198 //! _complex_float_view). Call .matrix to use as a GSL matrix (no copy is
199 //! involved). Allows one to use all GSL built-in functions. Note: non-owning
200 //! pointer - matrix AND gsl_view must remain in scope.
201 [[nodiscard]] auto as_gsl_view();
202
203 //! As above, but const
204 [[nodiscard]] auto as_gsl_view() const;
205
206 //============================================================================
207 // Basic matrix operations:
208 //============================================================================
209 //============================================================================
210
211 //! Returns the determinant. Uses GSL; via LU decomposition. Only works for
212 //! `double/complex<double>`
213 [[nodiscard]] T determinant() const;
214
215 //! Inverts the matrix, in place. Uses GSL; via LU decomposition. Only works
216 //! for `double/complex<double>.`
218
219 //! Returns inverse of the matrix. Leaves original matrix intact. Uses GSL;
220 //! via LU decomposition. Only works for `double/complex<double>`
221 [[nodiscard]] Matrix<T> inverse() const;
222 //! Returns transpose of matrix
223 [[nodiscard]] Matrix<T> transpose() const;
224
225 //============================================================================
226
227 //! Constructs a diagonal unit matrix (identity), in place; only for square
229 //! Sets all elements to zero, in place
230 Matrix<T> &zero();
231 // //! M -> M + aI, for I=identity (adds a to diag elements), in place
232 // Matrix<T> &plusIdent(T a = T(1));
233 //============================================================================
234
235 //! Returns conjugate of matrix
236 [[nodiscard]] Matrix<T> conj() const;
237 //! Returns real part of complex matrix (changes type; returns a real matrix)
238 [[nodiscard]] auto real() const;
239 //! Returns imag part of complex matrix (changes type; returns a real matrix)
240 [[nodiscard]] auto imag() const;
241 //! Converts a real to complex matrix (changes type; returns a complex matrix)
242 [[nodiscard]] auto complex() const;
243
244 //! Conjugates matrix, in place
246
247 //============================================================================
248 //! Muplitplies all the elements by those of matrix a, in place: M_ij *= a_ij
250
251 //! Returns new matrix, C_ij = A_ij*B_ij
252 [[nodiscard]] friend Matrix<T> mult_elements(Matrix<T> a,
253 const Matrix<T> &b) {
254 return a.mult_elements_by(b);
255 }
256
257 //============================================================================
258 // Operator overloads: +,-, scalar */
259 //! Overload standard operators: do what expected
260 Matrix<T> &operator+=(const Matrix<T> &rhs);
261 Matrix<T> &operator-=(const Matrix<T> &rhs);
262 Matrix<T> &operator*=(const T x);
263 Matrix<T> &operator/=(const T x);
264
265 //============================================================================
266 // nb: these are defined inline here to avoid ambiguous overload?
267 [[nodiscard]] friend Matrix<T> operator+(Matrix<T> lhs,
268 const Matrix<T> &rhs) {
269 return (lhs += rhs);
270 }
271 [[nodiscard]] friend Matrix<T> operator-(Matrix<T> lhs,
272 const Matrix<T> &rhs) {
273 return (lhs -= rhs);
274 }
275 [[nodiscard]] friend Matrix<T> operator*(const T x, Matrix<T> rhs) {
276 return (rhs *= x);
277 }
278 [[nodiscard]] friend Matrix<T> operator*(Matrix<T> lhs, const T x) {
279 return (lhs *= x);
280 }
281 [[nodiscard]] friend Matrix<T> operator/(Matrix<T> lhs, const T x) {
282 return (lhs /= x);
283 }
284
285 //============================================================================
286 //! Matrix<T> += T : T assumed to be *Identity!
287 Matrix<T> &operator+=(T aI);
288 //! Matrix<T> -= T : T assumed to be *Identity!
289 Matrix<T> &operator-=(T aI);
290
291 [[nodiscard]] friend Matrix<T> operator+(Matrix<T> M, T aI) {
292 return (M += aI);
293 }
294 [[nodiscard]] friend Matrix<T> operator-(Matrix<T> M, T aI) {
295 return (M -= aI);
296 }
297
298 //============================================================================
299 //! Matrix multiplication: C_ij = sum_k A_ik*B_kj
300 template <typename U>
301 friend Matrix<U> operator*(const Matrix<U> &a, const Matrix<U> &b);
302
303 //============================================================================
304 template <typename U>
305 friend std::ostream &operator<<(std::ostream &os, const Matrix<U> &a);
306};
307
308//==============================================================================
309//==============================================================================
310//==============================================================================
311
312//! Matrix multiplication C = A*B. C is overwritten with product, and must be correct size already
313/*! @details Wrapper for BLAS gemm functions:
314*/
315template <typename T>
316void GEMM(const Matrix<T> &a, const Matrix<T> &b, Matrix<T> *c,
317 bool trans_A = false, bool trans_B = false);
318
319//------------------------------------------------------------------------------
320
321//! 5-matrix contraction for N*N matrices: M_ab = A_ai B_aj C_ij D_ib E_jb, with BLAS
322/*!
323 @details
324 All input matrices are assumed square with the same dimension N.
325 The contraction sums over indices i and j:
326 \f[
327 M_{ab} = \sum_{i,j} A_{ai}\,B_{aj}\,C_{ij}\,D_{ib}\,E_{jb}.
328 \f]
329
330 This implementation evaluates the contraction by looping over the
331 shared index i and using a BLAS GEMM for the inner j-summation:
332 \f[
333 X^{(i)}_{jb} = C_{ij}\,E_{jb}, \qquad
334 Y^{(i)}_{ab} = \sum_j B_{aj}\,X^{(i)}_{jb},
335 \f]
336 followed by the accumulation
337 \f[
338 M_{ab} \mathrel{+}= A_{ai}\,Y^{(i)}_{ab}\,D_{ib}.
339 \f]
340
341 Thus each i-iteration performs one dense matrix multiplication
342 (B * X^{(i)}) via BLAS, while the remaining operations are
343 elementwise scalings and accumulations.
344
345 The routine allocates temporary work matrices X and Y of size NN
346 (reused across i) and writes the result into @p M, which must be
347 preallocated and is overwritten on entry (i.e., initialised to zero).
348
349 @note
350 ~20% faster than Naiive implementation, BUT not easily parallelisable.
351
352 @tparam T Scalar type (float, double, std::complex<...>)
353 @param A,B,C,D,E Input NN matrices
354 @param[out] M Output NN matrix receiving the contraction. Must be previously sized
355*/
356template <typename T>
357void PENTA_GEMM(const Matrix<T> &A, const Matrix<T> &B, const Matrix<T> &C,
358 const Matrix<T> &D, const Matrix<T> &E, Matrix<T> *M);
359
360//! M_ab = A_ai B_aj C_ij D_ib E_jb. Assuming all square matrices. Uses blas.
361/*!
362 @details
363 Calls:
364 \ref PENTA_GEMM(const Matrix<T>&, const Matrix<T>&,
365 const Matrix<T>&, const Matrix<T>&,
366 const Matrix<T>&, Matrix<T>*).
367*/
368template <typename T>
369Matrix<T> PENTA_GEMM(const Matrix<T> &A, const Matrix<T> &B, const Matrix<T> &C,
370 const Matrix<T> &D, const Matrix<T> &E) {
371 Matrix<T> M(A.rows(), A.cols());
372 PENTA_GEMM<T>(A, B, C, D, E, &M);
373 return M;
374}
375
376//-----------
377
378//! 5-matrix contraction for N*N matrices: M_ab = A_ai B_aj C_ij D_ib E_jb, without BLAS, but with optional OpenMP parallelisation
379/*! @details
380 Computes the tensor contraction
381 \f[
382 M_{ab} = \sum_{i,j} A_{ai}\,B_{aj}\,C_{ij}\,D_{ib}\,E_{jb},
383 \f]
384 assuming all input matrices are square with the same dimension.
385
386 This implementation uses explicit nested loops (no BLAS).
387 If @tparam PARALLEL is true, OpenMP pragmas are enabled to
388 parallelise the outer loops.
389
390 For a BLAS-accelerated implementation, see
391 \ref PENTA_GEMM(const Matrix<T>&, const Matrix<T>&,
392 const Matrix<T>&, const Matrix<T>&,
393 const Matrix<T>&, Matrix<T>*).
394
395 Usually, this version will be significantly faster than the BLAS one.
396 However, if it iself is being called in a parallel region, the BLAS one will be faster.
397
398 @tparam T Scalar type: double, float, std::complex
399 @tparam PARALLEL Enable OpenMP parallelisation when true (compile-time)
400 @param A,B,C,D,E Input NN matrices
401 @param[out] M Output NN matrix receiving the contraction. Must be previously sized
402*/
403template <typename T, bool PARALLEL = false>
404void PENTA(const Matrix<T> &A, const Matrix<T> &B, const Matrix<T> &C,
405 const Matrix<T> &D, const Matrix<T> &E, Matrix<T> *M);
406
407//! 5-index contraction for N*N matrices: M_ab = A_ai B_aj C_ij D_ib E_jb, without BLAS, but with optional OpenMP parallelisation
408//! @details calls above
409template <typename T, bool PARALLEL = false>
410Matrix<T> PENTA(const Matrix<T> &A, const Matrix<T> &B, const Matrix<T> &C,
411 const Matrix<T> &D, const Matrix<T> &E) {
412 Matrix<T> M(A.rows(), A.cols());
413 PENTA<T, PARALLEL>(A, B, C, D, E, &M);
414 return M;
415}
416
417//==============================================================================
418//==============================================================================
419//==============================================================================
420
421template <typename T>
422constexpr auto myEps();
423//! Compares two matrices; returns true iff all elements compare relatively to
424//! better than eps
425template <typename T>
426bool equal(const Matrix<T> &lhs, const Matrix<T> &rhs, T eps = myEps<T>());
427
428} // namespace LinAlg
429
430#include "Matrix.ipp"
Provides a non-owning view onto a Matrix (read/write, but no resize)
Definition Matrix.ipp:45
Matrix class; row-major.
Definition Matrix.hpp:39
T at(std::size_t row_i, std::size_t col_j) const
const () index access (with range checking). (i,j) ith row, jth col
Definition Matrix.hpp:140
const T * row(std::size_t row) const
Returns raw c pointer to start of a row.
Definition Matrix.hpp:166
Matrix< T > & make_identity()
Constructs a diagonal unit matrix (identity), in place; only for square.
Definition Matrix.ipp:201
Matrix< T > & invert_in_place()
Inverts the matrix, in place. Uses GSL; via LU decomposition. Only works for double/complex<double>.
Definition Matrix.ipp:135
const T & atc(std::size_t row_i, std::size_t col_j) const
const ref () index access (with range checking). (i,j) ith row, jth col
Definition Matrix.hpp:146
Matrix(std::size_t rows, std::size_t cols, const std::vector< T > &v)
Initialise a matrix from single initialiser list. {...}.
Definition Matrix.hpp:88
std::size_t size() const
Return rows*columns [total array size].
Definition Matrix.hpp:116
const T * data() const
As above, but const.
Definition Matrix.hpp:125
Matrix_view< const T > matrix_view() const
Returns an immutable view of the entire matrix (no ownership, no resize)
Definition Matrix.hpp:192
void resize(std::size_t rows, std::size_t cols, const T &value)
Resizes matrix to new dimension; all values reset to 'value'.
Definition Matrix.hpp:103
T * data()
Returns pointer to first element. Note: for std::complex<T>, this is a pointer to complex<T>,...
Definition Matrix.hpp:123
std::size_t rows() const
Return rows [major index size].
Definition Matrix.hpp:112
friend Matrix< U > operator*(const Matrix< U > &a, const Matrix< U > &b)
Matrix multiplication: C_ij = sum_k A_ik*B_kj.
Matrix< T > conj() const
Returns conjugate of matrix.
Definition Matrix.ipp:221
View< T > row_view(std::size_t row)
Returns a mutable 'View' of a row.
Definition Matrix.hpp:171
Matrix< T > & operator+=(const Matrix< T > &rhs)
Overload standard operators: do what expected.
Definition Matrix.ipp:276
auto begin()
iterators for underlying std::vector (entire data)
Definition Matrix.hpp:159
T & at(std::size_t row_i, std::size_t col_j)
() index access (with range checking). (i,j) returns ith row, jth col
Definition Matrix.hpp:135
auto complex() const
Converts a real to complex matrix (changes type; returns a complex matrix)
Definition Matrix.ipp:263
auto imag() const
Returns imag part of complex matrix (changes type; returns a real matrix)
Definition Matrix.ipp:252
auto as_gsl_view()
Returns gsl_matrix_view (or _float_view, _complex_view, _complex_float_view). Call ....
Definition Matrix.ipp:536
Matrix< T > & conj_in_place()
Conjugates matrix, in place.
Definition Matrix.ipp:232
T * operator[](std::size_t i)
As above, but const.
Definition Matrix.hpp:132
Matrix(std::size_t dimension)
Initialise a blank square matrix dimension*dimension, filled with 0.
Definition Matrix.hpp:60
Matrix()
Default initialiser.
Definition Matrix.hpp:48
Matrix(std::size_t rows, std::size_t cols, std::initializer_list< T > l)
Initialise a matrix from single initialiser list. {...}.
Definition Matrix.hpp:75
T & operator()(std::size_t i, std::size_t j)
() index access (with range checking). (i,j) returns ith row, jth col
Definition Matrix.hpp:152
Matrix< T > transpose() const
Returns transpose of matrix.
Definition Matrix.ipp:169
Matrix< T > & zero()
Sets all elements to zero, in place.
Definition Matrix.ipp:212
Matrix< T > inverse() const
Returns inverse of the matrix. Leaves original matrix intact. Uses GSL; via LU decomposition....
Definition Matrix.ipp:162
View< const T > row_view(std::size_t row) const
Returns an immutable 'View' of a row.
Definition Matrix.hpp:175
bool empty() const
Bool - is empty? (same as size==0)
Definition Matrix.hpp:119
std::size_t cols() const
Return columns [minor index size].
Definition Matrix.hpp:114
View< T > column_view(std::size_t col)
Returns a mutable 'View' of a column.
Definition Matrix.hpp:179
T determinant() const
Returns the determinant. Uses GSL; via LU decomposition. Only works for double/complex<double>
Definition Matrix.ipp:107
Matrix< T > & mult_elements_by(const Matrix< T > &a)
Muplitplies all the elements by those of matrix a, in place: M_ij *= a_ij.
Definition Matrix.ipp:328
const T * operator[](std::size_t i) const
[] index access (with no range checking). [i][j] returns ith row, jth col
Definition Matrix.hpp:130
auto real() const
Returns real part of complex matrix (changes type; returns a real matrix)
Definition Matrix.ipp:241
Matrix(std::initializer_list< std::initializer_list< T > > ll)
Initialise a matrix from initialiser list. {{},{},{}}. Each row must be same length.
Definition Matrix.hpp:65
Matrix_view< T > matrix_view()
Returns a mutable view of the entire matrix (no ownership, no resize)
Definition Matrix.hpp:188
Matrix(std::size_t rows, std::size_t cols, std::vector< T > &&v)
Initialise a matrix from single initialiser list. {...}.
Definition Matrix.hpp:82
T operator()(std::size_t i, std::size_t j) const
As above, but const.
Definition Matrix.hpp:154
Matrix(std::size_t rows, std::size_t cols)
Initialise a blank matrix rows*cols, filled with 0.
Definition Matrix.hpp:51
friend Matrix< T > mult_elements(Matrix< T > a, const Matrix< T > &b)
Returns new matrix, C_ij = A_ij*B_ij.
Definition Matrix.hpp:252
Matrix(std::size_t rows, std::size_t cols, const T &value)
Initialise a matrix rows*cols, filled with 'value'.
Definition Matrix.hpp:55
void resize(std::size_t rows, std::size_t cols)
Resizes matrix to new dimension; all values reset to default.
Definition Matrix.hpp:96
View< const T > column_view(std::size_t col) const
Returns an immutable 'View' of a column.
Definition Matrix.hpp:183
Proved a "view" onto an array.
Definition Matrix.ipp:7
constexpr bool is_complex_v
User-defined type-trait: Checks whether T is a std::complex type.
Definition AdamsMoulton.hpp:109
Defines Matrix, Vector classes, and linear some algebra functions.
Definition Matrix.hpp:26
void PENTA(const Matrix< T > &A, const Matrix< T > &B, const Matrix< T > &C, const Matrix< T > &D, const Matrix< T > &E, Matrix< T > *M)
5-matrix contraction for N*N matrices: M_ab = A_ai B_aj C_ij D_ib E_jb, without BLAS,...
Definition Matrix.ipp:488
void PENTA_GEMM(const Matrix< T > &A, const Matrix< T > &B, const Matrix< T > &C, const Matrix< T > &D, const Matrix< T > &E, Matrix< T > *M)
5-matrix contraction for N*N matrices: M_ab = A_ai B_aj C_ij D_ib E_jb, with BLAS
Definition Matrix.ipp:454
void GEMM(const Matrix< T > &a, const Matrix< T > &b, Matrix< T > *c, bool trans_A=false, bool trans_B=false)
Matrix multiplication C = A*B. C is overwritten with product, and must be correct size already.
Definition Matrix.ipp:394
bool equal(const Matrix< T > &lhs, const Matrix< T > &rhs, T eps=myEps< T >())
Compares two matrices; returns true iff all elements compare relatively to better than eps.
Definition Matrix.ipp:607