HILA
Loading...
Searching...
No Matches
matrix.h
Go to the documentation of this file.
1/**
2 * @file matrix.h
3 * @brief Definition of Matrix types
4 * @details This file contains base matrix type Matrix_t which defines all general matrix type
5 * operations Matrix types are Matrix, #Vector, #RowVector, #SquareMatrix of which Matrix is defined
6 * as a class and the rest are special cases of the Matrix class.
7 *
8 */
9
10#ifndef MATRIX_H_
11#define MATRIX_H_
12
13#include <type_traits>
14#include <sstream>
15#include "plumbing/defs.h"
16#include "datatypes/cmplx.h"
17
18// forward definitions of needed classes
19template <const int n, const int m, typename T, typename Mtype>
20class Matrix_t;
21
22template <const int n, const int m, typename T = double>
23class Array;
24
25template <int n, int m, typename T>
26class Matrix;
27
28template <int n, typename T>
29class DiagonalMatrix;
30
31/**
32 * @brief Vector is defined as 1-column Matrix
33 */
34template <int n, typename T>
36
37/**
38 * @brief RowVector is a 1-row Matrix
39 */
40template <int n, typename T>
42
43/**
44 * @brief Square matrix is defined as alias with special case of Matrix<n,n,T>
45 */
46template <int n, typename T>
48
49// template <const int n, const int m, typename T>
50// class DaggerMatrix;
51
52// Special case - m.column(), column of a matrix (not used now)
53// #include "matrix_column.h"
54
55
56/// @brief type to store the return combo of svd:
57/// {U, D, V} where U and V are nxn unitary / orthogonal,
58/// and D is real diagonal singular value matrices.
59/// @tparam M - type of input matrix
60template <typename M>
61struct svd_result {
62 static_assert(M::is_matrix() && M::rows() == M::columns(), "SVD only for square matrix");
63 M U;
64 DiagonalMatrix<M::size(), hila::arithmetic_type<M>> singularvalues;
65 M V;
66};
67
68/// @brief type to store the return value of eigen_hermitean():
69/// {E, U} where E is nxn DiagonalMatrix containing eigenvalues and
70/// U nxn unitary matrix, with eigenvector columns
71/// @tparam M - type of input matrix
72template <typename M>
74 static_assert(M::is_matrix() && M::rows() == M::columns(),
75 "Eigenvalues only for square matrix");
76 DiagonalMatrix<M::size(), hila::arithmetic_type<M>> eigenvalues;
77 M eigenvectors;
78};
79
80/**
81 * @brief The main \f$ n \times m \f$ matrix type template Matrix_t. This is a base class type for
82 * "useful" types which are derived from this.
83 *
84 * @details Uses curiously recurring template pattern (CRTP), where the last template parameter is
85 * the template itself
86 *
87 * Example: the Matrix type below is defined as
88 * @code{.cpp}
89 * template <int n, int m, typename T>
90 * class Matrix : public Matrix_t<n, m, T, Matrix<n, m, T>> { .. }
91 * @endcode
92 *
93 * This pattern is used because stupid c++ makes it complicated to write generic code, in this case
94 * derived functions to return derived type
95 *
96 * @tparam n Number of rows
97 * @tparam m Number of columns
98 * @tparam T Matrix element type
99 * @tparam Mtype Specific "Matrix" type for CRTP
100 */
101
102// Helper struct for getting the floating point number epsilons without
103// having to use std::numeric_limits .
104
105template <const int n, const int m, typename T, typename Mtype>
106class Matrix_t {
107
108 public:
109 /// The data as a one dimensional array
110 T c[n * m];
111
112 public:
114 "Matrix requires Complex or arithmetic type");
115
116 // std incantation for field types
117 using base_type = hila::arithmetic_type<T>;
118 using argument_type = T;
119
120 // help for templates, can use T::is_matrix()
121 // Not very useful outside template parameters
122 static constexpr bool is_matrix() {
123 return true;
124 }
125
126 /**
127 * @brief Returns true if Matrix is a vector
128 *
129 * @return true
130 * @return false
131 */
132 static constexpr bool is_vector() {
133 return (n == 1 || m == 1);
134 }
135
136 /**
137 * @brief Returns true if matrix is a square matrix
138 *
139 * @return true
140 * @return false
141 */
142 static constexpr bool is_square() {
143 return (n == m);
144 }
145
146 /// Define default constructors to ensure std::is_trivial
147 // Move constructor not needed
148 Matrix_t() = default;
149 ~Matrix_t() = default;
150 Matrix_t(const Matrix_t &v) = default;
151
152 // constructor from scalar -- keep it explicit! Not good for auto use
153 // NOTE: I forgot why this should be kept explicit
154 template <typename S, int nn = n, int mm = m,
155 std::enable_if_t<(hila::is_assignable<T &, S>::value && nn == mm), int> = 0>
156 explicit inline Matrix_t(const S rhs) {
157 for (int i = 0; i < n; i++)
158 for (int j = 0; j < n; j++) {
159 if (i == j)
160 e(i, j) = rhs;
161 else
162 e(i, j) = 0;
163 }
164 }
165
166 // Construct from a different type matrix
167 template <typename S, typename MT,
168 std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
169 Matrix_t(const Matrix_t<n, m, S, MT> &rhs) out_only {
170 for (int i = 0; i < n * m; i++) {
171 c[i] = rhs.c[i];
172 }
173 }
174
175 // construct from 0
176 inline Matrix_t(const std::nullptr_t &z) {
177 for (int i = 0; i < n * m; i++) {
178 c[i] = 0;
179 }
180 }
181
182 // Construct matrix automatically from right-size initializer list
183 // This does not seem to be dangerous, so keep non-explicit
184 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
185 inline Matrix_t(std::initializer_list<S> rhs) {
186 assert(rhs.size() == n * m &&
187 "Matrix/Vector initializer list size must match variable size");
188 int i = 0;
189 for (auto it = rhs.begin(); it != rhs.end(); it++, i++) {
190 c[i] = *it;
191 }
192 }
193
194 // cast to curious type
195
196 template <typename Tm = Mtype,
197 std::enable_if_t<!std::is_same<Tm, Matrix<n, m, T>>::value, int> = 0>
198 inline operator Mtype &() {
199 return *reinterpret_cast<Mtype *>(this);
200 }
201
202
203 template <typename Tm = Mtype,
204 std::enable_if_t<!std::is_same<Tm, Matrix<n, m, T>>::value, int> = 0>
205 inline operator const Mtype &() const {
206 return *reinterpret_cast<const Mtype *>(this);
207 }
208
209 // automatically cast to generic matrix
210
211 inline operator Matrix<n, m, T> &() {
212 return *reinterpret_cast<Matrix<n, m, T> *>(this);
213 }
214
215 inline operator const Matrix<n, m, T> &() const {
216 return *reinterpret_cast<const Matrix<n, m, T> *>(this);
217 }
218
219 /// Define constant methods rows(), columns() - may be useful in template code
220
221 /**
222 * @brief Returns row length
223 *
224 * @return constexpr int
225 */
226 static constexpr int rows() {
227 return n;
228 }
229 /**
230 * @brief Returns column length
231 *
232 * @return constexpr int
233 */
234 static constexpr int columns() {
235 return m;
236 }
237
238
239 /**
240 * @brief Returns size of #Vector or square Matrix
241 *
242 * @tparam q row size n
243 * @tparam p column size m
244 * @return constexpr int
245 */
246 // size for row vector
247 template <int q = n, int p = m, std::enable_if_t<q == 1, int> = 0>
248 static constexpr int size() {
249 return p;
250 }
251 // size for column vector
252 template <int q = n, int p = m, std::enable_if_t<p == 1, int> = 0>
253 static constexpr int size() {
254 return q;
255 }
256 // size for square matrix
257 template <int q = n, int p = m, std::enable_if_t<q == p, int> = 0>
258 static constexpr int size() {
259 return q;
260 }
261
262 /**
263 * @brief Standard array indexing operation for matrices and vectors
264 *
265 * @details Accessing singular elements is insufficient, but matrix elements are often quite
266 * small.
267 *
268 * Exammple for matrix:
269 * \code
270 * Matrix<n,m,MyType> M;
271 * MyType a = M.e(i,j); \\ i <= n, j <= m
272 * \endcode
273 *
274 * Example for vector:
275 * \code {.cpp}
276 * Vector<n,MyType> V;
277 * MyType a = V.e(i) \\ i <= n
278 * \endcode
279 *
280 * @param i row index
281 * @param j column index
282 * @return T matrix element type
283 */
284
285#pragma hila loop_function
286 inline T e(const int i, const int j) const {
287 // return elem[i][j];
288 return c[i * m + j];
289 }
290 // Same as above but with const_function, see const_function for details
291 inline T &e(const int i, const int j) const_function {
292 // return elem[i][j];
293 return c[i * m + j];
294 }
295 // declare single e here too in case we have a vector
296 // (n || m == 1)
297
298#pragma hila loop_function
299 template <int q = n, int p = m, std::enable_if_t<(q == 1 || p == 1), int> = 0>
300 inline T e(const int i) const {
301 return c[i];
302 }
303 // Same as above but with const_function, see const_function for details
304 template <int q = n, int p = m, std::enable_if_t<(q == 1 || p == 1), int> = 0>
305 inline T &e(const int i) const_function {
306 return c[i];
307 }
308
309 /**
310 * @brief Indexing operation [] defined only for vectors.
311 *
312 * @details Example:
313 *
314 * \code {.cpp}
315 * Vector<n,MyType> V;
316 * MyType a = V[i] \\ i <= n
317 * \endcode
318 *
319 * @tparam q row size n
320 * @tparam p column size m
321 * @param i row or vector index depending on which is being indexed
322 * @return T
323 */
324
325#pragma hila loop_function
326 template <int q = n, int p = m, std::enable_if_t<(q == 1 || p == 1), int> = 0>
327 inline T operator[](const int i) const {
328 return c[i];
329 }
330 // Same as above but with const_function, see const_function for details
331 template <int q = n, int p = m, std::enable_if_t<(q == 1 || p == 1), int> = 0>
332 inline T &operator[](const int i) const_function {
333 return c[i];
334 }
335
336 /**
337 * @brief Return row of a matrix as a RowVector
338 *
339 * @param r index of row to be referenced
340 * @return const RowVector<m, T>&
341 */
342 RowVector<m, T> row(int r) const {
344 for (int i = 0; i < m; i++)
345 v[i] = e(r, i);
346 return v;
347 }
348
349 /**
350 * @brief Set row of Matrix with #RowVector if types are assignable
351 *
352 * @tparam S RowVector type
353 * @param r Index of row to be set
354 * @param v RowVector to be set
355 */
356 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
357 void set_row(int r, const RowVector<m, S> &v) {
358 for (int i = 0; i < m; i++)
359 e(r, i) = v[i];
360 }
361
362 /**
363 * @brief Returns column vector as value at index c
364 *
365 * @param c index of column vector to be returned
366 * @return const Vector<n, T>
367 */
368 Vector<n, T> column(int c) const {
369 Vector<n, T> v;
370 for (int i = 0; i < n; i++)
371 v[i] = e(i, c);
372 return v;
373 }
374
375
376 // matrix_col_t<n,T,Mtype> column(int c) {
377
378 // return matrix_col_t<n,T,Mtype>(*this,c);
379 // }
380
381
382 /// get column of a matrix
383 // hila_matrix_column_t<n, T, Mtype> column(int c) {
384 // return hila_matrix_column_t<n, T, Mtype>(*this, c);
385 // }
386
387 /**
388 * @brief Set column of Matrix with #Vector if types are assignable
389 *
390 * @tparam S Vector type
391 * @param c Index of column to be set
392 * @param v #Vector to be set
393 */
394 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
395 void set_column(int c, const Vector<n, S> &v) {
396 for (int i = 0; i < n; i++)
397 e(i, c) = v[i];
398 }
399
400 /**
401 * @brief Return diagonal of square matrix
402 * @details If called for non square matrix the program will throw an error.
403 *
404 * @return Vector<n, T> returned vector.
405 */
407 static_assert(n == m, "diagonal() method defined only for square matrices");
409 for (int i = 0; i < n; i++)
410 res.e(i) = (*this).e(i, i);
411 return res;
412 }
413
414 /**
415 * @brief Set diagonal of square matrix to #Vector which is passed to the method
416 * @details If called for non square matrix the program will throw an error.
417 *
418 * Example:
419 *
420 * \code {.cpp}
421 * SquareMatrix<n,MyType> S = 0; \\ Zero matrix
422 * Vector<n,MyType> V = 1; \\ Vector assigned to 1 at all elements
423 * S.set_diagonal(V); \\ Results in Identity matrix of size n
424 * \endcode
425 *
426 * @tparam S type vector to assign values to
427 * @param v Vector to assign to diagonal
428 */
429 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
430 void set_diagonal(const Vector<n, S> &v) {
431 static_assert(n == m, "set_diagonal() method defined only for square matrices");
432 for (int i = 0; i < n; i++)
433 (*this).e(i, i) = v.e(i);
434 }
435
436
437 /**
438 * @brief Cast Matrix to Array
439 * @details used for array operations
440 *
441 * @return Array<n, m, T>&
442 */
443 const Array<n, m, T> &asArray() const {
444 return *reinterpret_cast<const Array<n, m, T> *>(this);
445 }
446 // Same as above but with const_function, see const_function for details
447 Array<n, m, T> &asArray() const_function {
448 return *reinterpret_cast<Array<n, m, T> *>(this);
449 }
450
451 /**
452 * @brief Cast Vector to DiagonalMatrix
453 *
454 * @return DiagonalMatrix<n,T>
455 */
456 #pragma hila loop_function
457 template <int mm = m, std::enable_if_t<mm == 1, int> = 0>
459 return *reinterpret_cast<const DiagonalMatrix<n, T> *>(this);
460 }
461 // Same as above but with const_function, see const_function for details
462 template <int mm = m, std::enable_if_t<mm == 1, int> = 0>
463 DiagonalMatrix<n, T> &asDiagonalMatrix() const_function {
464 return *reinterpret_cast<DiagonalMatrix<n, T> *>(this);
465 }
466
467
468 /// casting from one Matrix (number) type to another: do not do this automatically.
469 /// but require an explicit cast operator. This makes it easier to write code.
470 /// or should it be automatic? keep/remove explicit?
471 /// TODO: CHECK AVX CONVERSIONS
472
473 // template <typename S, typename Rtype,
474 // std::enable_if_t<
475 // Rtype::is_matrix() && Rtype::rows() == n && Rtype::columns() == m
476 // &&
477 // hila::is_assignable<typename (Rtype::argument_type) &,
478 // T>::value,
479 // int> = 0>
480 // explicit operator Rtype() const {
481 // Rtype res;
482 // for (int i = 0; i < n * m; i++)
483 // res.c[i] = c[i];
484 // return res;
485 // }
486
487 /**
488 * @brief Unary - operator
489 * @details Returns matrix with the signs of all the elements in the Matrix flipped.
490 *
491 * @return Mtype
492 */
493 inline Mtype operator-() const {
494 Mtype res;
495 for (int i = 0; i < n * m; i++) {
496 res.c[i] = -c[i];
497 }
498 return res;
499 }
500
501 /**
502 * @brief Unary + operator
503 * @details Equivalent to identity operator meaning that matrix stays as is.
504 *
505 * @return const Mtype&
506 */
507 inline const auto &operator+() const {
508 return *this;
509 }
510
511 /**
512 * @brief Boolean operator == to determine if two matrices are exactly the same
513 * @details Tolerance for equivalence is zero, meaning that the matrices must be exactly the
514 * same.
515 *
516 * @tparam S Type for Matrix which is being compared to
517 * @param rhs right hand side Matrix which we are comparing
518 * @return true
519 * @return false
520 */
521 template <typename S, int n2, int m2>
522 bool operator==(const Matrix<n2, m2, S> &rhs) const {
523 if constexpr (n != n2 || m != m2)
524 return false;
525
526 for (int i = 0; i < n; i++)
527 for (int j = 0; j < m; j++) {
528 if (e(i, j) != rhs.e(i, j))
529 return false;
530 }
531 return true;
532 }
533
534 /**
535 * @brief Boolean operator != to check if matrices are exactly different
536 * @details if matrices are exactly the same then this will return false
537 *
538 * @tparam S Type for MAtrix which is being compared to
539 * @param rhs right hand side Matrix which we are comparing
540 * @return true
541 * @return false
542 */
543 template <typename S>
544 bool operator!=(const Matrix<n, m, S> &rhs) const {
545 return !(*this == rhs);
546 }
547
548 /**
549 * @brief Assignment operator = to assign values to matrix
550 * @details The following ways to assign a matrix are:
551 *
552 *
553 * __Assignment from Matrix__:
554 *
555 * \code {.cpp}
556 * Matrix<n,m,MyType> M_0;
557 * .
558 * . M_0 has values assigned to it
559 * .
560 * Matrix<n,m,MyType> M; \\ undefined matrix
561 * M = M_0; \\ Assignment from M_0
562 * \endcode
563 *
564 * __Assignment from 0__:
565 *
566 * \code {.cpp}
567 * Matrix<n,m,MyType> M;
568 * M = 0; Zero matrix;
569 * \endcode
570 *
571 * __Assignment from scalar__:
572 *
573 * Assignment from scalar assigns the scalar to the diagonal elements as \f$ M = I\cdot a\f$
574 *
575 * \code {.cpp}
576 * MyType a = hila::random;
577 * Matrix<n,m,MyType> M;
578 * M = a; M = I*a
579 * \endcode
580 *
581 *__Initializer list__:
582 *
583 * Assignment from c++ initializer list.
584 *
585 * \code{.cpp}
586 * Matrix<2,2,int> M ;
587 * M = {1, 0
588 * 0, 1};
589 * \endcode
590 */
591
592 template <typename S, typename MT,
593 std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
594 inline auto &operator=(const Matrix_t<n, m, S, MT> &rhs) out_only & {
595 for (int i = 0; i < n * m; i++) {
596 c[i] = rhs.c[i];
597 }
598 return *this;
599 }
600
601 // assign from 0
602
603 inline auto &operator=(const std::nullptr_t &z) out_only & {
604 for (int i = 0; i < n * m; i++) {
605 c[i] = 0;
606 }
607 return *this;
608 }
609
610 // Assign from "scalar" for square matrix
611
612 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value && n == m, int> = 0>
613 inline auto &operator=(const S rhs) out_only & {
614
615 for (int i = 0; i < n; i++)
616 for (int j = 0; j < n; j++) {
617 if (i == j)
618 e(i, j) = rhs;
619 else
620 e(i, j) = 0;
621 }
622 return *this;
623 }
624
625 // Assign from diagonal matrix
626
627 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
628 inline auto &operator=(const DiagonalMatrix<n, S> &rhs) out_only & {
629 static_assert(n == m,
630 "Assigning DiagonalMatrix to Matrix possible only for square matrices");
631
632 for (int i = 0; i < n; i++)
633 for (int j = 0; j < n; j++) {
634 if (i == j)
635 e(i, j) = rhs.e(i);
636 else
637 e(i, j) = 0;
638 }
639 return *this;
640 }
641
642
643 // Assign from initializer list
644
645 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
646 auto &operator=(std::initializer_list<S> rhs) out_only & {
647 assert(rhs.size() == n * m && "Initializer list has a wrong size in assignment");
648 int i = 0;
649 for (auto it = rhs.begin(); it != rhs.end(); it++, i++) {
650 c[i] = *it;
651 }
652 return *this;
653 }
654
655 // Delete the rvalue-assign op
656 template <typename S>
657 Matrix_t &operator=(const S &s) && = delete;
658
659
660 /**
661 * @brief Add assign operator with matrix or scalar
662 * @details Add assign operator can be used in the following ways
663 *
664 * __Add assign matrix__:
665 *
666 * \code {.cpp}
667 * Matrix<n,m,MyType> M,N;
668 * M = 1;
669 * N = 1;
670 * M += N; \\M = 2*I
671 * \endcode
672 *
673 * __Add assign scalar__:
674 *
675 * Adds scalar \f$ a \f$ to __square__ matrix as \f$ M + a\cdot\mathbb{1} \f$
676 *
677 * \code {.cpp}
678 * Matrix<n,m,MyType> M = 1;
679 * M += 1 ; \\ M = 2*I
680 * \endcode
681 *
682 * @tparam S Element type of rhs
683 * @tparam MT Matrix type of rhs
684 * @param rhs Matrix to add
685 * @return Mtype&
686 */
687
688 template <typename S, typename MT,
689 std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
690 auto &operator+=(const Matrix_t<n, m, S, MT> &rhs) & {
691 for (int i = 0; i < n * m; i++) {
692 c[i] += rhs.c[i];
693 }
694 return *this;
695 }
696
697 /**
698 * @brief Subtract assign operator with matrix or scalar
699 * @details Subtract assign operator can be used in the following ways
700 *
701 * __Subtract assign matrix__:
702 *
703 * \code {.cpp}
704 * Matrix<n,m,MyType> M,N;
705 * M = 3;
706 * N = 1;
707 * M -= N; \\M = 2*I
708 * \endcode
709 *
710 * __Subtract assign scalar__:
711 *
712 * Subtract scalar \f$ a \f$ to __square__ matrix as \f$ M - a\cdot\mathbb{1} \f$
713 *
714 * \code {.cpp}
715 * Matrix<n,m,MyType> M = 3;
716 * M -= 1 ; \\ M = 2*I
717 * \endcode
718 *
719 * @param rhs Matrix to subtract with
720 * @return template <typename S, typename MT,
721 * std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>&
722 */
723
724 template <typename S, typename MT,
725 std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
726 auto &operator-=(const Matrix_t<n, m, S, MT> &rhs) & {
727 for (int i = 0; i < n * m; i++) {
728 c[i] -= rhs.c[i];
729 }
730 return *this;
731 }
732
733 // add assign a scalar to square matrix
734
735 template <typename S,
736 std::enable_if_t<hila::is_assignable<T &, hila::type_plus<T, S>>::value, int> = 0>
737 auto &operator+=(const S &rhs) & {
738
739 static_assert(n == m, "rows != columns : scalar addition possible for square matrix only!");
740
741 for (int i = 0; i < n; i++) {
742 e(i, i) += rhs;
743 }
744 return *this;
745 }
746
747 // subtract assign type T and convertible
748
749 template <typename S,
750 std::enable_if_t<hila::is_assignable<T &, hila::type_minus<T, S>>::value, int> = 0>
751 auto &operator-=(const S rhs) & {
752 static_assert(n == m,
753 "rows != columns : scalar subtraction possible for square matrix only!");
754 for (int i = 0; i < n; i++) {
755 e(i, i) -= rhs;
756 }
757 return *this;
758 }
759
760 /**
761 * @brief Multiply assign scalar or matrix
762 * @details Multiplication works as defined for matrix multiplication and scalar matrix
763 * multiplication.
764 *
765 * Matrix multiply assign only defined for square matrices, since the matrix dimensions would
766 * change otherwise.
767 *
768 * Multiply assign operator can be used in the following ways
769 *
770 * __Multiply assign matrix__:
771 *
772 * \code {.cpp}
773 * Matrix<n,m,MyType> M,N;
774 * .
775 * . Fill matrices M and N
776 * .
777 * M *= N; \\ M = M*N
778 * \endcode
779 *
780 * __Multiply assign scalar__:
781 *
782 * \code {.cpp}
783 * Matrix<n,m,MyType> M;
784 * .
785 * . Fill whole matrix with 1
786 * .
787 * M *= 2 ; \\ M is filled with 2
788 * \endcode
789 *
790 * @param rhs Matrix to multiply with
791 * @return template <int p, typename S, typename MT,
792 * std::enable_if_t<hila::is_assignable<T &, hila::type_mul<T, S>>::value, int> = 0>&
793 */
794 template <int p, typename S, typename MT,
795 std::enable_if_t<hila::is_assignable<T &, hila::type_mul<T, S>>::value, int> = 0>
796 auto &operator*=(const Matrix_t<m, p, S, MT> &rhs) & {
797 static_assert(m == p, "can't assign result of *= to lhs Matrix, because doing so "
798 "would change it's dimensions");
799 *this = *this * rhs;
800 return *this;
801 }
802
803
804 /*
805 // same type square matrices:
806 template <int p,typename S,typename MT,
807 std::enable_if_t<hila::is_assignable<T&,hila::type_mul<T,S>>::value, int> = 0>
808 Mtype& operator*=(const Matrix_t<m,p,S,MT>& rhs)& {
809 static_assert(m==p,"can't assign result of *= to lhs Matrix, because doing so "
810 "would change it's dimensions");
811
812 S tmp_row[m];
813 int i,j,k;
814 for(i=0; i<m; ++i) {
815 for(j=0; j<m; ++j) {
816 tmp_row[j]=e(i,j);
817 }
818 for(j=0; j<m; ++j) {
819 e(i,j)=tmp_row[0]*rhs.e(0,j);
820 for(k=1; k<m; ++k) {
821 e(i,j)+=tmp_row[k]*rhs.e(k,j);
822 }
823 }
824 }
825 return *this;
826 }
827 */
828
829 // multiply assign with scalar
830
831 template <typename S,
832 std::enable_if_t<hila::is_assignable<T &, hila::type_mul<T, S>>::value, int> = 0>
833 auto &operator*=(const S rhs) & {
834 for (int i = 0; i < n * m; i++) {
835 c[i] *= rhs;
836 }
837 return *this;
838 }
839
840 /**
841 * @brief Divide assign scalar
842 * @details Divide works as defined for scalar matrix division.
843 *
844 * Division assign operator can be used in the following ways
845 *
846 * __Divide assign scalar__:
847 *
848 * \code {.cpp}
849 * Matrix<n,m,MyType> M;
850 * .
851 * . Fill whole matrix with 2
852 * .
853 * M /= 2 ; \\ M is filled with 1
854 * \endcode
855 *
856 * @param rhs Matrix to divide with
857 * @return template <int p, typename S, typename MT,
858 * std::enable_if_t<hila::is_assignable<T &, hila::type_mul<T, S>>::value, int> = 0>&
859 */
860
861 template <typename S,
862 std::enable_if_t<hila::is_assignable<T &, hila::type_div<T, S>>::value, int> = 0>
863 auto &operator/=(const S rhs) & {
864 for (int i = 0; i < n * m; i++) {
865 c[i] /= rhs;
866 }
867 return *this;
868 }
869
870 /**
871 * @brief add and sub assign a DiagonalMatrix
872 *
873 * This is possible only for square matrices
874 */
875 template <typename S,
876 std::enable_if_t<hila::is_assignable<T &, hila::type_plus<T, S>>::value, int> = 0>
877 auto &operator+=(const DiagonalMatrix<n, S> &rhs) & {
878 static_assert(n == m, "Assigning DiagonalMatrix possible only for square matrix");
879
880 for (int i = 0; i < n; i++)
881 e(i, i) += rhs.e(i);
882 return *this;
883 }
884
885 template <typename S,
886 std::enable_if_t<hila::is_assignable<T &, hila::type_plus<T, S>>::value, int> = 0>
887 auto &operator-=(const DiagonalMatrix<n, S> &rhs) & {
888 static_assert(n == m, "Assigning DiagonalMatrix possible only for square matrix");
889
890 for (int i = 0; i < n; i++)
891 e(i, i) -= rhs.e(i);
892 return *this;
893 }
894
895 /**
896 * @brief mult and divide assign a diagonal - cols must match diagonal matrix rows
897 */
898 template <typename S,
899 std::enable_if_t<hila::is_assignable<T &, hila::type_mul<T, S>>::value, int> = 0>
900 auto &operator*=(const DiagonalMatrix<m, S> &rhs) & {
901
902 for (int i = 0; i < n; i++)
903 for (int j = 0; j < m; j++)
904 e(i, j) *= rhs.e(j);
905
906 return *this;
907 }
908
909 template <typename S,
910 std::enable_if_t<hila::is_assignable<T &, hila::type_div<T, S>>::value, int> = 0>
911 auto &operator/=(const DiagonalMatrix<m, S> &rhs) & {
912
913 for (int i = 0; i < n; i++)
914 for (int j = 0; j < m; j++)
915 e(i, j) /= rhs.e(j);
916
917 return *this;
918 }
919
920
921 /**
922 * @brief Matrix fill
923 * @details Fills the matrix with element if it is assignable to matrix type T
924 *
925 * Works as follows:
926 *
927 * \code {.cpp}
928 * Matrix<n,m,MyType> M;
929 * M.fill(2) \\ Matrix is filled with 2
930 * \endcode
931 *
932 * @tparam S Element type to be assigned
933 * @param rhs Element to fill matrix with
934 * @return const Mtype&
935 */
936 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
937 const auto &fill(const S rhs) out_only {
938 for (int i = 0; i < n * m; i++)
939 c[i] = rhs;
940 return *this;
941 }
942
943 /**
944 * @brief Transpose of matrix
945 * @details Return type for square matrix is same input, for non square return type is
946 * Matrix<n,m,MyType>
947 * @tparam mm
948 * @tparam std::conditional<n, Mtype, Matrix<m, n, T>>::type
949 * @return const Rtype
950 */
951 template <int mm = m,
952 typename Rtype = typename std::conditional<n == m, Matrix_t, Matrix<m, n, T>>::type,
953 std::enable_if_t<(mm != 1 && n != 1), int> = 0>
954 inline Rtype transpose() const {
955 Rtype res;
956 for (int i = 0; i < n; i++)
957 for (int j = 0; j < m; j++) {
958 res.e(j, i) = e(i, j);
959 }
960 return res;
961 }
962
963 /**
964 * @brief Transpose of vector
965 * @details Returns reference
966 *
967 * @tparam mm
968 * @return const RowVector<n, T>&
969 */
970 template <int mm = m, std::enable_if_t<mm == 1, int> = 0>
971 inline const RowVector<n, T> &transpose() const {
972 return *reinterpret_cast<const RowVector<n, T> *>(this);
973 }
974
975 template <int nn = n, std::enable_if_t<nn == 1 && m != 1, int> = 0>
976 inline const Vector<m, T> &transpose() const {
977 return *reinterpret_cast<const Vector<m, T> *>(this);
978 }
979
980
981 /**
982 * @brief Returns complex conjugate of Matrix
983 *
984 * @return const Mtype
985 */
986 inline Mtype conj() const {
987 Mtype res;
988 for (int i = 0; i < n * m; i++) {
989 res.c[i] = ::conj(c[i]);
990 }
991 return res;
992 }
993
994 /**
995 * @brief Hermitian conjugate of matrix
996 * @details for square matrix return type is same, for non square it is Matrix<m,n,MyType>
997 *
998 * @tparam std::conditional<n, Mtype, Matrix<m, n, T>>::type
999 * @return const Rtype
1000 */
1001 template <typename Rtype = typename std::conditional<n == m, Mtype, Matrix<m, n, T>>::type>
1002 inline Rtype dagger() const {
1003 Rtype res;
1004 for (int i = 0; i < n; i++)
1005 for (int j = 0; j < m; j++) {
1006 res.e(j, i) = ::conj(e(i, j));
1007 }
1008 return res;
1009 }
1010
1011 /**
1012 * @brief Adjoint of matrix
1013 * @details Alias to dagger
1014 *
1015 * @tparam std::conditional<n, Mtype, Matrix<m, n, T>>::type
1016 * @return Rtype
1017 */
1018 template <typename Rtype = typename std::conditional<n == m, Mtype, Matrix<m, n, T>>::type>
1019 inline Rtype adjoint() const {
1020 return dagger();
1021 }
1022
1023 /**
1024 * @brief Returns absolute value of Matrix
1025 * @details For Matrix<n,m,Complex<T>> case type is changed to Matrix<n,m,T> as expected since
1026 * absolute value of a complex number is real
1027 *
1028 * @tparam M
1029 * @return Mtype
1030 */
1031 auto abs() const {
1033 for (int i = 0; i < n * m; i++) {
1034 res.c[i] = ::abs(c[i]);
1035 }
1036 return res;
1037 }
1038
1039
1040 // It seems that using special "Dagger" type makes the code slower!
1041 // Disable it now
1042 // inline const DaggerMatrix<m,n,T> & dagger() const {
1043 // return *reinterpret_cast<const DaggerMatrix<m,n,T>*>(this);
1044 // }
1045
1046 /**
1047 * @brief Returns real part of Matrix or #Vector
1048 *
1049 * @return Matrix<n, m, hila::arithmetic_type<T>>
1050 */
1053 for (int i = 0; i < m * n; i++) {
1054 res.c[i] = ::real(c[i]);
1055 }
1056 return res;
1057 }
1058
1059 /**
1060 * @brief Returns imaginary part of Matrix or #Vector
1061 *
1062 * @return Matrix<n, m, hila::arithmetic_type<T>>
1063 */
1066 for (int i = 0; i < m * n; i++) {
1067 res.c[i] = ::imag(c[i]);
1068 }
1069 return res;
1070 }
1071
1072 /**
1073 * @brief Computes Trace for Matrix
1074 * @details Not define for non square matrices. Static assert will stop code from compiling if
1075 * executed for non square matrices.
1076 *
1077 * @return T
1078 */
1079 T trace() const {
1080 static_assert(n == m, "trace not defined for non square matrices");
1081 T result(0);
1082 for (int i = 0; i < n; i++) {
1083 result += e(i, i);
1084 }
1085 return result;
1086 }
1087
1088 /**
1089 * @brief Multiply with given matrix and compute trace of result
1090 * @details Slightly cheaper operation than \code (M*N).trace() \endcode
1091 *
1092 * @tparam p
1093 * @tparam q
1094 * @tparam S
1095 * @tparam MT
1096 * @param rm
1097 * @return hila::type_mul<T, S>
1098 */
1099 template <int p, int q, typename S, typename MT>
1100 hila::type_mul<T, S> mul_trace(const Matrix_t<p, q, S, MT> &rm) const {
1101
1102 static_assert(p == m && q == n, "mul_trace(): argument matrix size mismatch");
1103
1104 hila::type_mul<T, S> res = 0;
1105 for (int i = 0; i < n; i++)
1106 for (int j = 0; j < m; j++) {
1107 res += e(i, j) * rm.e(j, i);
1108 }
1109 return res;
1110 }
1111
1112 /**
1113 * @brief Calculate square norm - sum of squared elements
1114 *
1115 * @return hila::arithmetic_type<T>
1116 */
1117 hila::arithmetic_type<T> squarenorm() const {
1118 hila::arithmetic_type<T> result(0);
1119 for (int i = 0; i < n * m; i++) {
1120 result += ::squarenorm(c[i]);
1121 }
1122 return result;
1123 }
1124
1125 /**
1126 * @brief Calculate vector norm - sqrt of squarenorm
1127 *
1128 * @tparam S
1129 * @return hila::arithmetic_type<T>
1130 */
1131 template <typename S = T,
1132 std::enable_if_t<hila::is_floating_point<hila::arithmetic_type<S>>::value, int> = 0>
1133 hila::arithmetic_type<T> norm() const {
1134 return sqrt(squarenorm());
1135 }
1136
1137 template <typename S = T,
1138 std::enable_if_t<!hila::is_floating_point<hila::arithmetic_type<S>>::value, int> = 0>
1139 double norm() const {
1140 return sqrt(static_cast<double>(squarenorm()));
1141 }
1142
1143 /**
1144 * @brief Find max or min value - only for arithmetic types
1145 */
1146
1147 template <typename S = T, std::enable_if_t<hila::is_arithmetic<S>::value, int> = 0>
1148 T max() const {
1149 T res = c[0];
1150 for (int i = 1; i < n * m; i++) {
1151 if (res < c[i])
1152 res = c[i];
1153 }
1154 return res;
1155 }
1156
1157 template <typename S = T, std::enable_if_t<hila::is_arithmetic<S>::value, int> = 0>
1158 T min() const {
1159 T res = c[0];
1160 for (int i = 1; i < n * m; i++) {
1161 if (res > c[i])
1162 res = c[i];
1163 }
1164 return res;
1165 }
1166
1167 auto max_abs() const {
1168 hila::arithmetic_type<T> tres, res = 0;
1169 for (int i = 0; i < n * m; i++) {
1170 tres = ::abs(c[i]);
1171 if (tres > res) {
1172 res = tres;
1173 }
1174 }
1175 return res;
1176 }
1177
1178 auto min_abs() const {
1179 hila::arithmetic_type<T> tres, res = ::abs(c[0]);
1180 for (int i = 1; i < n * m; i++) {
1181 tres = ::abs(c[i]);
1182 if (tres < res) {
1183 res = tres;
1184 }
1185 }
1186 return res;
1187 }
1188
1189 // dot product - (*this).dagger() * rhs
1190 // could be done as well by writing the operation as above!
1191 /**
1192 * @brief Dot product
1193 * @details Only works between two #Vector objects
1194 *
1195 * \code {.cpp}
1196 * Vector<m,MyType> V,W;
1197 * .
1198 * .
1199 * .
1200 * V.dot(W); // valid operation
1201 * \endcode
1202 *
1203 * \code {.cpp}
1204 * RowVector<m,MyType> V,W;
1205 * .
1206 * .
1207 * .
1208 * V.dot(W); // not valid operation
1209 * \endcode
1210 *
1211 *
1212 * @tparam p Row length for rhs
1213 * @tparam q Column length for rhs
1214 * @tparam S Type for rhs
1215 * @tparam R Gives resulting type of lhs and rhs multiplication
1216 * @param rhs Vector to compute dot product with
1217 * @return R Value of dot product
1218 */
1219 template <int p, int q, typename S, typename R = hila::type_mul<T, S>>
1220 inline R dot(const Matrix<p, q, S> &rhs) const {
1221 static_assert(m == 1 && q == 1 && p == n,
1222 "dot() product only for vectors of the same length");
1223
1224 R r = 0;
1225 for (int i = 0; i < n; i++) {
1226 r += ::conj(c[i]) * rhs.e(i);
1227 }
1228 return r;
1229 }
1230
1231 // outer product - (*this) * rhs.dagger(), sizes (n,1) and (p,1)
1232 // gives n * p matrix
1233 // could be done as well by the above operation!
1234
1235 /**
1236 * @brief Outer product
1237 * @details Only works between two #Vector objects
1238 *
1239 * \code{.cpp}
1240 * Vector<n,MyType> V,W;
1241 * .
1242 * .
1243 * .
1244 * V.outer_product(W); \\ Valid operation
1245 * \endcode
1246 *
1247 * \code {.cpp}
1248 * RowVector<m,MyType> V,W;
1249 * .
1250 * .
1251 * .
1252 * V.outer_product(W); // not valid operation
1253 * \endcode
1254 *
1255 * @tparam p Row length for rhs
1256 * @tparam q Column length for rhs
1257 * @tparam S Element type for rhs
1258 * @tparam R Type between lhs and rhs multiplication
1259 * @param rhs Vector to compute outer product with
1260 * @return Matrix<n, p, R>
1261 */
1262 template <int p, int q, typename S, typename R = hila::type_mul<T, S>>
1264 static_assert(m == 1 && q == 1, "outer_product() only for vectors");
1265
1266 Matrix<n, p, R> res;
1267 for (int i = 0; i < n; i++) {
1268 for (int j = 0; j < p; j++) {
1269 res.e(i, j) = c[i] * ::conj(rhs.e(j));
1270 }
1271 }
1272 return res;
1273 }
1274
1275
1276 /// dot with matrix - matrix
1277 // template <int p, std::enable_if_t<(p > 1 || m > 1), int> = 0>
1278 // inline Matrix<m, p, T> dot(const Matrix<n, p, T> &rhs) const {
1279 // Matrix<m, p, T> res;
1280 // for (int i = 0; i < m; i++)
1281 // for (int j = 0; j < p; j++) {
1282 // res.e(i, j) = 0;
1283 // for (int k = 0; k < n; j++)
1284 // res.e(i, j) += ::conj(e(k, i)) * rhs.e(k, j);
1285 // }
1286 // }
1287
1288 /// Generate random elements
1289
1290 /**
1291 * @brief Fills Matrix with random elements
1292 * @details Works only for real valued elements such as float or double
1293 *
1294 * @return Mtype&
1295 */
1296 Mtype &random() out_only {
1297
1298 static_assert(hila::is_floating_point<hila::arithmetic_type<T>>::value,
1299 "Matrix/Vector random() requires non-integral type elements");
1300
1301 for (int i = 0; i < n * m; i++) {
1302 hila::random(c[i]);
1303 }
1304 return *this;
1305 }
1306
1307 /**
1308 * @brief Fills Matrix with gaussian random elements
1309 * @details Works only for real valued elements such as float or double
1310 *
1311 * @param width
1312 * @return Mtype&
1313 */
1314 Mtype &gaussian_random(base_type width = 1.0) out_only {
1315
1316 static_assert(hila::is_floating_point<hila::arithmetic_type<T>>::value,
1317 "Matrix/Vector gaussian_random() requires non-integral type elements");
1318
1319 // for Complex numbers gaussian_random fills re and im efficiently
1320 if constexpr (hila::is_complex<T>::value) {
1321 for (int i = 0; i < n * m; i++) {
1322 c[i].gaussian_random(width);
1323 }
1324 } else {
1325 // now not complex matrix
1326 // if n*m even, max i in loop below is n*m-2.
1327 // if n*m odd, max i is n*m-3
1328 double gr;
1329 for (int i = 0; i < n * m - 1; i += 2) {
1330 c[i] = hila::gaussrand2(gr) * width;
1331 c[i + 1] = gr * width;
1332 }
1333 if constexpr ((n * m) % 2 > 0) {
1334 c[n * m - 1] = hila::gaussrand() * width;
1335 }
1336 }
1337 return *this;
1338 }
1339
1340 /**
1341 * @brief permute columns of Matrix
1342 * @details Reordering is done based off of permutation vector
1343 *
1344 * @param permutation Vector of integers to permute columns
1345 * @return Mtype
1346 */
1347 Mtype permute_columns(const Vector<m, int> &permutation) const {
1348 Mtype res;
1349 for (int i = 0; i < m; i++)
1350 res.set_column(i, this->column(permutation[i]));
1351 return res;
1352 }
1353
1354 /**
1355 * @brief permute rows of Matrix
1356 * @details Reordering is done based off of permutation vector
1357 *
1358 * @param permutation Vector of integers to permute rows
1359 * @return Mtype
1360 */
1361 Mtype permute_rows(const Vector<n, int> &permutation) const {
1362 Mtype res;
1363 for (int i = 0; i < n; i++)
1364 res.set_row(i, this->row(permutation[i]));
1365 return res;
1366 }
1367
1368
1369 /**
1370 * @brief Permute elements of vector
1371 * @details Reordering is done based off of permutation vector
1372 *
1373 * @param permutation Vector of integers to permute rows
1374 * @return Mtype
1375 */
1376 template <int N>
1377 Mtype permute(const Vector<N, int> &permutation) const {
1378 static_assert(
1379 n == 1 || m == 1,
1380 "permute() only for vectors, use permute_rows() or permute_columns() for matrices");
1381 static_assert(N == Mtype::size(), "Incorrect size of permutation vector");
1382
1383 Mtype res;
1384 for (int i = 0; i < N; i++) {
1385 res[i] = (*this)[permutation[i]];
1386 }
1387 return res;
1388 }
1389
1390 /**
1391 * @brief Sort method for #Vector
1392 * @details Two interfaces: first returns permutation vector, which can be used to permute
1393 * other vectors/matrices second does only sort
1394 *
1395 * __Direct sort__:
1396 *
1397 * @code {.cpp}
1398 * Vector<n,MyType> V;
1399 * V.random();
1400 * V.sort(); // V is sorted
1401 * @endcode
1402 *
1403 * __Permutation vector__:
1404 *
1405 * @code {.cpp}
1406 * Vector<n,MyType> V;
1407 * Vector<n,int> perm;
1408 * V.random();
1409 * V.sort(perm);
1410 * V.permute(perm); // V is sorted
1411 * @endcode
1412 *
1413 */
1414#pragma hila novector
1415 template <int N>
1416 Mtype sort(Vector<N, int> &permutation, hila::sort order = hila::sort::ascending) const {
1417
1418 static_assert(n == 1 || m == 1, "Sorting possible only for vectors");
1419 static_assert(hila::is_arithmetic<T>::value,
1420 "Sorting possible only for arithmetic vector elements");
1421 static_assert(N == Mtype::size(), "Incorrect size of permutation vector");
1422
1423 for (int i = 0; i < N; i++)
1424 permutation[i] = i;
1425 if (hila::sort::unsorted == order) {
1426 return *this;
1427 }
1428
1429 if (hila::sort::ascending == order) {
1430 for (int i = 1; i < N; i++) {
1431 for (int j = i; j > 0 && c[permutation[j - 1]] > c[permutation[j]]; j--)
1432 hila::swap(permutation[j], permutation[j - 1]);
1433 }
1434 } else {
1435 for (int i = 1; i < N; i++) {
1436 for (int j = i; j > 0 && c[permutation[j - 1]] < c[permutation[j]]; j--)
1437 hila::swap(permutation[j], permutation[j - 1]);
1438 }
1439 }
1440
1441 return this->permute(permutation);
1442 }
1443
1444#pragma hila novector
1445 Mtype sort(hila::sort order = hila::sort::ascending) const {
1446 static_assert(n == 1 || m == 1, "Sorting possible only for vectors");
1447
1448 Vector<Mtype::size(), int> permutation;
1449 return sort(permutation, order);
1450 }
1451
1452
1453 /**
1454 * @brief Multiply \f$ n \times m \f$-matrix from the left by \f$ n \times m \f$ matrix defined
1455 * by \f$ 2 \times 2 \f$ sub matrix
1456 * @details The multiplication is defined as follows, let \f$M\f$ as the \f$ 2 \times 2 \f$
1457 * input matrix and \f$B\f$ be `(this)` matrix, being the matrix stored in the object this
1458 * method is called for. Let \f$A = I\f$ be a \f$ n \times m \f$ unit matrix. We then set the
1459 * values of A to be: \f{align}{ A_{p,p} = M_{0,0}, \hspace{5px} A_{p,q} = M_{0,1}, \hspace{5px}
1460 * A_{q,p} = M_{1,0}, \hspace{5px} A_{q,q} = M_{1,1}. \f}
1461 *
1462 * Then the resulting matrix will be:
1463 *
1464 * \f{align}{ B = A \cdot B \f}
1465 *
1466 * @tparam R Element type of M
1467 * @tparam Mt Matrix type of M
1468 * @param p First row and column
1469 * @param q Second row and column
1470 * @param M \f$ 2 \times 2\f$ Matrix to multiply with
1471 */
1472 template <typename R, typename Mt>
1473 void mult_by_2x2_left(int p, int q, const Matrix_t<2, 2, R, Mt> &M) {
1474 static_assert(hila::contains_complex<T>::value || !hila::contains_complex<R>::value,
1475 "Cannot multiply real matrix with complex 2x2 block matrix");
1476
1477 Vector<2, T> a;
1478 for (int i = 0; i < m; i++) {
1479 a.e(0) = this->e(p, i);
1480 a.e(1) = this->e(q, i);
1481
1482 a = M * a;
1483
1484 this->e(p, i) = a.e(0);
1485 this->e(q, i) = a.e(1);
1486 }
1487 }
1488
1489 /**
1490 * @brief Multiply \f$ n \times m \f$-matrix from the right by \f$ n \times m \f$ matrix
1491 * defined by \f$ 2 \times 2 \f$ sub matrix
1492 * @details See Matrix::mult_by_2x2_left, only difference being that the multiplication is from
1493 * the right.
1494 *
1495 * @tparam R Element type of M
1496 * @tparam Mt Matrix type of M
1497 * @param p First row and column
1498 * @param q Second row and column
1499 * @param M \f$ 2 \times 2\f$ Matrix to multiply with
1500 */
1501 template <typename R, typename Mt>
1502 void mult_by_2x2_right(int p, int q, const Matrix_t<2, 2, R, Mt> &M) {
1503 static_assert(hila::contains_complex<T>::value || !hila::contains_complex<R>::value,
1504 "Cannot multiply real matrix with complex 2x2 block matrix");
1505
1507 for (int i = 0; i < n; i++) {
1508 a.e(0) = this->e(i, p);
1509 a.e(1) = this->e(i, q);
1510
1511 a = a * M;
1512
1513 this->e(i, p) = a.e(0);
1514 this->e(i, q) = a.e(1);
1515 }
1516 }
1517
1518
1519 // Following methods are defined in matrix_linalg.h
1520 //
1521
1522 /**
1523 * @brief following calculate the determinant of a square matrix
1524 * det() is the generic interface, using laplace for small matrices and LU for large
1525 */
1526
1527#pragma hila novector
1528 T det_lu() const;
1529#pragma hila novector
1530 T det_laplace() const;
1531#pragma hila novector
1532 T det() const;
1533
1534 /**
1535 * @brief Calculate eigenvalues and -vectors of an hermitean (or symmetric) matrix.
1536 * @details Returns the number of Jacobi iterations, or -1 if did not converge.
1537 * Arguments will contain eigenvalues and eigenvectors in columns of the "eigenvectors" matrix.
1538 * Computation is done in double precision despite the input matrix types.
1539 * @param eigenvaluevec Vector of computed eigenvalue
1540 * @param eigenvectors Eigenvectors as columns of \f$ n \times n \f$ Matrix
1541 *
1542 */
1543
1544#pragma hila novector
1545 template <typename Et, typename Mt, typename MT>
1546 int eigen_hermitean(out_only DiagonalMatrix<n, Et> &eigenvalues,
1547 out_only Matrix_t<n, n, Mt, MT> &eigenvectors,
1548 enum hila::sort sorted = hila::sort::unsorted) const;
1549
1550 eigen_result<Mtype> eigen_hermitean(enum hila::sort sorted = hila::sort::unsorted) const;
1551
1552 // Temporary interface to eigenvalues, to be removed
1553 template <typename Et, typename Mt, typename MT>
1554 int eigen_jacobi(out_only Vector<n, Et> &eigenvaluevec,
1555 out_only Matrix_t<n, n, Mt, MT> &eigenvectors,
1556 enum hila::sort sorted = hila::sort::unsorted) const {
1557
1559 int i = eigen_hermitean(d, eigenvectors, sorted);
1560 eigenvaluevec = d.asArray().asVector();
1561 return i;
1562 }
1563
1564
1565#pragma hila novector
1566 template <typename Et, typename Mt, typename MT>
1567 int svd(out_only Matrix_t<n, n, Mt, MT> &_U, out_only DiagonalMatrix<n, Et> &_D,
1568 out_only Matrix_t<n, n, Mt, MT> &_V,
1569 enum hila::sort sorted = hila::sort::unsorted) const;
1570
1571 svd_result<Mtype> svd(enum hila::sort sorted = hila::sort::unsorted) const;
1572
1573
1574#pragma hila novector
1575 template <typename Et, typename Mt, typename MT>
1576 int svd_pivot(out_only Matrix_t<n, n, Mt, MT> &_U, out_only DiagonalMatrix<n, Et> &_D,
1577 out_only Matrix_t<n, n, Mt, MT> &_V,
1578 enum hila::sort sorted = hila::sort::unsorted) const;
1579
1580 svd_result<Mtype> svd_pivot(enum hila::sort sorted = hila::sort::unsorted) const;
1581
1582
1583 //////// matrix_linalg.h
1584
1585 /// Convert to string for printing
1586 ///
1587 std::string str(int prec = 8, char separator = ' ') const {
1588 std::stringstream text;
1589 text.precision(prec);
1590 for (int i = 0; i < n; i++) {
1591 for (int j = 0; j < m; j++) {
1592 text << e(i, j);
1593 if (i < n - 1 || j < m - 1)
1594 text << separator;
1595 }
1596 }
1597 return text.str();
1598 }
1599};
1600
1601/**
1602 * @brief \f$ n \times m \f$ Matrix class which defines matrix operations.
1603 * @details The Matrix class is a derived class of Matrix_t, which is the general definition of a
1604 * matrix class. See Matrix_t details section for reasoning.
1605 *
1606 * All mathematical operations for Matrix are inherited from Matrix_t and are visible on this page.
1607 * __NOTE__: Some method documentation is not being displayed, for example the assignment operator
1608 * documentation is not being inherited. See Matrix_t::operator= for use.
1609 *
1610 * To see all methods of initializing a matrix see constructor method #Matrix::Matrix
1611 *
1612 * __NOTE__: In the documentation examples n,m are integers and MyType is a HILA [standard
1613 * type](@ref standard) or Complex.
1614 *
1615 * @tparam n Number of rows
1616 * @tparam m Number of columns
1617 * @tparam T Data type Matrix
1618 */
1619template <int n, int m, typename T>
1620class Matrix : public Matrix_t<n, m, T, Matrix<n, m, T>> {
1621
1622 public:
1623 /// std incantation for field types
1624 using base_type = hila::arithmetic_type<T>;
1625 using argument_type = T;
1626
1627 /**
1628 * @brief Construct a new Matrix object
1629 * @details The following ways of constructing a matrix are:
1630 *
1631 * __NOTE__: n,m are integers and MyType is a HILA [standard type](@ref standard) or Complex.
1632 *
1633 * __Default constructor__:
1634 *
1635 * Allocates undefined \f$ n\times m\f$ Array.
1636 *
1637 * \code{.cpp}
1638 * Matrix<n,m,MyType> M;
1639 * \endcode
1640 *
1641 *
1642 * __Scalar constructor__:
1643 *
1644 * Construct with given scalar at diagonal elements \f$ M = \mathbf{I}\cdot x\f$.
1645 *
1646 * \code{.cpp}
1647 * MyType x = hila::random();
1648 * Matrix<n,m,MyType> M = x;
1649 * \endcode
1650 *
1651 * __Copy constructor__:
1652 *
1653 * Construction from previously defined matrix if types are compatible. For example the code
1654 * below only works if assignment from MyOtherType to MyType is defined.
1655 *
1656 * \code{.cpp}
1657 * Matrix<n,m,MyOtherType> M_0;
1658 * //
1659 * // M_0 is filled with content
1660 * //
1661 * Matrix<n,m,MyType> M = M_0;
1662 * \endcode
1663 *
1664 * __Initializer list__:
1665 *
1666 * Construction from c++ initializer list.
1667 *
1668 * \code{.cpp}
1669 * Matrix<2,2,int> M = {1, 0
1670 * 0, 1};
1671 * \endcode
1672 */
1673
1674 // use the Base::Base -trick to inherit constructors and assignments
1675 using Matrix_t<n, m, T, Matrix<n, m, T>>::Matrix_t;
1676
1677 using Matrix_t<n, m, T, Matrix<n, m, T>>::operator=;
1678 Matrix &operator=(const Matrix &v) out_only & = default;
1679};
1680
1681namespace hila {
1682
1683//////////////////////////////////////////////////////////////////////////
1684// Tool to get "right" result type for matrix (T1) + (T2) -op, where
1685// T1 and T2 are either complex or arithmetic matrices
1686// Use: hila::mat_x_mat_type<T1,T2>
1687// - T1 == T2 returns T1
1688// - If T1 or T2 contains complex, returns Matrix<rows,cols,Complex<typeof T1+T2>>
1689// - otherwise returns the "larger" type
1690
1691template <typename T1, typename T2>
1692struct matrix_op_type_s {
1693 using type = Matrix<T1::rows(), T1::columns(), hila::ntype_op<T1, T2>>;
1694};
1695
1696// if types are the same
1697template <typename T>
1698struct matrix_op_type_s<T, T> {
1699 using type = T;
1700};
1701
1702// useful format mat_x_mat_type<T1,T2>
1703template <typename T1, typename T2>
1704using mat_x_mat_type = typename matrix_op_type_s<T1, T2>::type;
1705
1706////////////////////////////////////////////////////////////////////////////////
1707// Matrix + scalar result type:
1708// hila::mat_scalar_type<Mt,S>
1709// - if result is convertible to Mt, return Mt
1710// - if Mt is not complex and S is, return
1711// Matrix<Complex<type_sum(scalar_type(Mt),scalar_type(S))>>
1712// - otherwise return Matrix<type_sum>
1713
1714template <typename Mt, typename S, typename Enable = void>
1715struct matrix_scalar_op_s {
1716 using type =
1717 Matrix<Mt::rows(), Mt::columns(),
1718 Complex<hila::type_plus<hila::arithmetic_type<Mt>, hila::arithmetic_type<S>>>>;
1719};
1720
1721template <typename Mt, typename S>
1722struct matrix_scalar_op_s<
1723 Mt, S,
1724 typename std::enable_if_t<std::is_convertible<hila::type_plus<hila::number_type<Mt>, S>,
1725 hila::number_type<Mt>>::value>> {
1726 // using type = Mt;
1727 using type = typename std::conditional<
1728 hila::is_floating_point<hila::arithmetic_type<Mt>>::value, Mt,
1729 Matrix<Mt::rows(), Mt::columns(),
1730 hila::type_plus<hila::arithmetic_type<Mt>, hila::arithmetic_type<S>>>>::type;
1731};
1732
1733template <typename Mt, typename S>
1734using mat_scalar_type = typename matrix_scalar_op_s<Mt, S>::type;
1735
1736
1737} // namespace hila
1738
1739//////////////////////////////////////////////////////////////////////////
1740
1741/**
1742 * @brief Return transposed Matrix or #Vector
1743 * @details Wrapper around Matrix::transpose
1744 *
1745 * Can be called as:
1746 * \code {.cpp}
1747 * Matrix<n,m,MyType> M,M_transposed;
1748 * M_transposed = transpose(M);
1749 * \endcode
1750 *
1751 *
1752 * @tparam Mtype Matrix type
1753 * @param arg Object to be transposed
1754 * @return auto Mtype transposed matrix
1755 */
1756template <typename Mtype, std::enable_if_t<Mtype::is_matrix(), int> = 0>
1757inline auto transpose(const Mtype &arg) {
1758 return arg.transpose();
1759}
1760
1761/**
1762 * @brief Return conjugate Matrix or #Vector
1763 * @details Wrapper around Matrix::conj
1764 *
1765 * Can be called as:
1766 * \code {.cpp}
1767 * Matrix<n,m,MyType> M,M_conjugate;
1768 * M_conjugate = conj(M);
1769 * \endcode
1770 *
1771 *
1772 * @tparam Mtype Matrix type
1773 * @param arg Object to be conjugated
1774 * @return auto Mtype conjugated matrix
1775 */
1776template <typename Mtype, std::enable_if_t<Mtype::is_matrix(), int> = 0>
1777inline auto conj(const Mtype &arg) {
1778 return arg.conj();
1779}
1780
1781/**
1782 * @brief Return adjoint Matrix
1783 * @details Wrapper around Matrix::adjoint
1784 *
1785 * Can be called as:
1786 * \code {.cpp}
1787 * Matrix<n,m,MyType> M,M_adjoint;
1788 * M_conjugate = adjoint(M);
1789 * \endcode
1790 *
1791 *
1792 * @tparam Mtype Matrix type
1793 * @param arg Object to compute adjoint of
1794 * @return auto Mtype adjoint matrix
1795 */
1796template <typename Mtype, std::enable_if_t<Mtype::is_matrix(), int> = 0>
1797inline auto adjoint(const Mtype &arg) {
1798 return arg.adjoint();
1799}
1800
1801/**
1802 * @brief Return dagger of Matrix
1803 * @details Wrapper around Matrix::adjoint
1804 *
1805 * Same as adjoint
1806 *
1807 *
1808 * @tparam Mtype Matrix type
1809 * @param arg Object to compute dagger of
1810 * @return auto Mtype dagger matrix
1811 */
1812template <typename Mtype, std::enable_if_t<Mtype::is_matrix(), int> = 0>
1813inline auto dagger(const Mtype &arg) {
1814 return arg.adjoint();
1815}
1816
1817/**
1818 * @brief Return absolute value Matrix or #Vector
1819 * @details Wrapper around Matrix::abs
1820 *
1821 * Can be called as:
1822 * \code {.cpp}
1823 * Matrix<n,m,MyType> M,M_abs;
1824 * M_abs = abs(M);
1825 * \endcode
1826 *
1827 *
1828 * @tparam Mtype Matrix type
1829 * @param arg Object to compute absolute value of
1830 * @return auto Mtype absolute value Matrix or Vector
1831 */
1832template <typename Mtype, std::enable_if_t<Mtype::is_matrix(), int> = 0>
1833inline auto abs(const Mtype &arg) {
1834 return arg.abs();
1835}
1836
1837/**
1838 * @brief Return trace of square Matrix
1839 * @details Wrapper around Matrix::trace
1840 *
1841 * Can be called as:
1842 * \code {.cpp}
1843 * Matrix<n,m,MyType> M;
1844 * MyType T;
1845 * T = trace(M);
1846 * \endcode
1847 *
1848 *
1849 * @tparam Mtype Matrix type
1850 * @param arg Object to compute trace of
1851 * @return auto Trace of Matrix
1852 */
1853template <typename Mtype, std::enable_if_t<Mtype::is_matrix(), int> = 0>
1854inline auto trace(const Mtype &arg) {
1855 return arg.trace();
1856}
1857
1858/**
1859 * @brief Return real of Matrix or #Vector
1860 * @details Wrapper around Matrix::real
1861 *
1862 * Can be called as:
1863 * \code {.cpp}
1864 * Matrix<n,m,MyType> M,M_real;
1865 * M_real = real(M);
1866 * \endcode
1867 *
1868 *
1869 * @tparam Mtype Matrix type
1870 * @param arg Object to compute real part of
1871 * @return auto Mtype real part of arg
1872 */
1873template <typename Mtype, std::enable_if_t<Mtype::is_matrix(), int> = 0>
1874inline auto real(const Mtype &arg) {
1875 return arg.real();
1876}
1877
1878/**
1879 * @brief Return imaginary of Matrix or #Vector
1880 * @details Wrapper around Matrix::imag
1881 *
1882 * Can be called as:
1883 * \code {.cpp}
1884 * Matrix<n,m,MyType> M,M_imag;
1885 * M_imag = imag(M);
1886 * \endcode
1887 *
1888 *
1889 * @tparam Mtype Matrix type
1890 * @param arg Object to compute imag part of
1891 * @return auto Mtype imag part of arg
1892 */
1893template <typename Mtype, std::enable_if_t<Mtype::is_matrix(), int> = 0>
1894inline auto imag(const Mtype &arg) {
1895 return arg.imag();
1896}
1897
1898
1899///////////////////////////////////////////////////////////////////////////
1900// Now matrix additions: matrix + matrix
1901
1902/**
1903 * @brief Addition operator
1904 * @details Defined for the following operations
1905 *
1906 * __Matrix + Matrix:__
1907 *
1908 * Addition operator between matrices is defined in the usual way (element wise).
1909 *
1910 * __NOTE__: Matrices must share same dimensionality.
1911 *
1912 * \code {.cpp}
1913 * Matrix<n,m,MyType> M, N, S;
1914 * M.fill(1);
1915 * N.fill(1);
1916 * S = M + N; // Resulting matrix is uniformly 2
1917 * \endcode
1918 *
1919 *
1920 * __Scalar + Matrix / Matrix + Scalar:__
1921 *
1922 * Addition operator between matrix and scalar is defined in the usual way, where the scalar is
1923 * added to the diagonal elements.
1924 *
1925 * __NOTE__: Exact definition exist in overloaded functions that can be viewed in source code.
1926 *
1927 * \f$ M + 2 = M + 2\cdot\mathbb{1} \f$
1928 *
1929 * \code {.cpp}
1930 * Matrix<n,m,MyType> M,S;
1931 * M.fill(0);
1932 * S = M + 1; // Resulting matrix is identity matrix
1933 * \endcode
1934 *
1935 *
1936 * @tparam Mtype1 Matrix type for a
1937 * @tparam Mtype2 Matrix type for b
1938 * @param a Left matrix or scalar
1939 * @param b Right matrix or scalar
1940 * @return Rtype Return matrix of compatible type between Mtype1 and Mtype2
1941 */
1942template <typename Mtype1, typename Mtype2,
1943 std::enable_if_t<Mtype1::is_matrix() && Mtype2::is_matrix(), int> = 0,
1944 typename Rtype = hila::mat_x_mat_type<Mtype1, Mtype2>,
1945 std::enable_if_t<!std::is_same<Mtype1, Rtype>::value, int> = 0>
1946inline Rtype operator+(const Mtype1 &a, const Mtype2 &b) {
1947
1948 constexpr int n = Mtype1::rows();
1949 constexpr int m = Mtype1::columns();
1950
1951 static_assert(n == Mtype2::rows() && m == Mtype2::columns(), "Matrix sizes do not match");
1952
1953 Rtype r;
1954 for (int i = 0; i < n; i++)
1955 for (int j = 0; j < m; j++)
1956 r.e(i, j) = a.e(i, j) + b.e(i, j);
1957 return r;
1958}
1959
1960// Real micro-optimization wrt. above - no extra creation of variable and copy.
1961
1962template <typename Mtype1, typename Mtype2,
1963 std::enable_if_t<Mtype1::is_matrix() && Mtype2::is_matrix(), int> = 0,
1964 typename Rtype = hila::mat_x_mat_type<Mtype1, Mtype2>,
1965 std::enable_if_t<std::is_same<Mtype1, Rtype>::value, int> = 0>
1966inline Rtype operator+(Mtype1 a, const Mtype2 &b) {
1967
1968 constexpr int n = Mtype1::rows();
1969 constexpr int m = Mtype1::columns();
1970
1971 static_assert(n == Mtype2::rows() && m == Mtype2::columns(), "Matrix sizes do not match");
1972
1973 for (int i = 0; i < n; i++)
1974 for (int j = 0; j < m; j++)
1975 a.e(i, j) += b.e(i, j);
1976 return a;
1977}
1978
1979/**
1980 * @brief Subtraction operator
1981 * @details Defined for the following operations
1982 *
1983 * __Matrix - Matrix:__
1984 *
1985 * Subtraction operator between matrices is defined in the usual way (element wise).
1986 *
1987 * __NOTE__: Matrices must share same dimensionality.
1988 *
1989 * \code {.cpp}
1990 * Matrix<n,m,MyType> M, N, S;
1991 * M.fill(1);
1992 * N.fill(1);
1993 * S = M - N; // Resulting matrix is uniformly 0
1994 * \endcode
1995 *
1996 *
1997 * __Scalar - Matrix / Matrix - Scalar:__
1998 *
1999 * Subtraction operator between matrix and scalar is defined in the usual way, where the scalar is
2000 * subtracted from the diagonal elements.
2001 *
2002 * __NOTE__: Exact definition exist in overloaded functions that can be viewed in source code.
2003 *
2004 * \f$ M - 2 = M - 2\cdot\mathbb{1} \f$
2005 *
2006 * \code {.cpp}
2007 * Matrix<n,m,MyType> M,S;
2008 * M = 2; // M = 2*I
2009 * S = M - 1; // Resulting matrix is identity matrix
2010 * \endcode
2011 *
2012 *
2013 * @tparam Mtype1 Matrix type for a
2014 * @tparam Mtype2 Matrix type for b
2015 * @param a Left matrix or scalar
2016 * @param b Right matrix or scalar
2017 * @return Rtype Return matrix of compatible type between Mtype1 and Mtype2
2018 */
2019template <typename Mtype1, typename Mtype2,
2020 std::enable_if_t<Mtype1::is_matrix() && Mtype2::is_matrix(), int> = 0,
2021 typename Rtype = hila::mat_x_mat_type<Mtype1, Mtype2>>
2022inline Rtype operator-(const Mtype1 &a, const Mtype2 &b) {
2023
2024 constexpr int n = Mtype1::rows();
2025 constexpr int m = Mtype1::columns();
2026
2027 static_assert(n == Mtype2::rows() && m == Mtype2::columns(), "Matrix sizes do not match");
2028
2029 Rtype r;
2030 for (int i = 0; i < n; i++)
2031 for (int j = 0; j < m; j++)
2032 r.e(i, j) = a.e(i, j) - b.e(i, j);
2033 return r;
2034}
2035
2036// Matrix + scalar
2037template <typename Mtype, typename S,
2038 std::enable_if_t<Mtype::is_matrix() && hila::is_complex_or_arithmetic<S>::value, int> = 0,
2039 typename Rtype = hila::mat_scalar_type<Mtype, S>>
2040inline Rtype operator+(const Mtype &a, const S &b) {
2041
2042 static_assert(Mtype::rows() == Mtype::columns(),
2043 "Matrix + scalar possible only for square matrix");
2044
2045 Rtype r;
2046 for (int i = 0; i < Rtype::rows(); i++)
2047 for (int j = 0; j < Rtype::columns(); j++) {
2048 r.e(i, j) = a.e(i, j);
2049 if (i == j)
2050 r.e(i, j) += b;
2051 }
2052 return r;
2053}
2054
2055// scalar + matrix
2056template <typename Mtype, typename S,
2057 std::enable_if_t<Mtype::is_matrix() && hila::is_complex_or_arithmetic<S>::value, int> = 0,
2058 typename Rtype = hila::mat_scalar_type<Mtype, S>>
2059inline Rtype operator+(const S &b, const Mtype &a) {
2060
2061 static_assert(Mtype::rows() == Mtype::columns(),
2062 "Matrix + scalar possible only for square matrix");
2063 return a + b;
2064}
2065
2066// matrix - scalar
2067template <typename Mtype, typename S,
2068 std::enable_if_t<Mtype::is_matrix() && hila::is_complex_or_arithmetic<S>::value, int> = 0,
2069 typename Rtype = hila::mat_scalar_type<Mtype, S>>
2070inline Rtype operator-(const Mtype &a, const S &b) {
2071
2072 static_assert(Mtype::rows() == Mtype::columns(),
2073 "Matrix - scalar possible only for square matrix");
2074
2075 Rtype r;
2076 for (int i = 0; i < Rtype::rows(); i++)
2077 for (int j = 0; j < Rtype::columns(); j++) {
2078 r.e(i, j) = a.e(i, j);
2079 if (i == j)
2080 r.e(i, j) -= b;
2081 }
2082 return r;
2083}
2084
2085// scalar - matrix
2086template <typename Mtype, typename S,
2087 std::enable_if_t<Mtype::is_matrix() && hila::is_complex_or_arithmetic<S>::value, int> = 0,
2088 typename Rtype = hila::mat_scalar_type<Mtype, S>>
2089inline Rtype operator-(const S &b, const Mtype &a) {
2090
2091 static_assert(Mtype::rows() == Mtype::columns(),
2092 "Scalar - matrix possible only for square matrix");
2093
2094 Rtype r;
2095 for (int i = 0; i < Rtype::rows(); i++)
2096 for (int j = 0; j < Rtype::columns(); j++) {
2097 r.e(i, j) = -a.e(i, j);
2098 if (i == j)
2099 r.e(i, j) += b;
2100 }
2101 return r;
2102}
2103
2104////////////////////////////////////////
2105// matrix * matrix is the most important bit
2106
2107// same type square matrices:
2108template <typename Mt, std::enable_if_t<Mt::is_matrix() && Mt::is_square(), int> = 0>
2109inline Mt operator*(const Mt &a, const Mt &b) {
2110
2111 constexpr int n = Mt::rows();
2112
2113 Mt res;
2114
2115 for (int i = 0; i < n; i++)
2116 for (int j = 0; j < n; j++) {
2117 res.e(i, j) = a.e(i, 0) * b.e(0, j);
2118 for (int k = 1; k < n; k++) {
2119 res.e(i, j) += a.e(i, k) * b.e(k, j);
2120 }
2121 }
2122 return res;
2123}
2124
2125/**
2126 * @brief Multiplication operator
2127 * @details Multiplication type depends on the original types of the multiplied matrices. Defined
2128 * for the following operations.
2129 *
2130 * __Matrix * Matrix:__
2131 *
2132 * [Standard matrix multiplication](https://en.wikipedia.org/wiki/Matrix_multiplication) where LHS
2133 * columns must match RHS rows
2134 *
2135 * \code {.cpp}
2136 * Matrix<2, 3, double> M;
2137 * Matrix<3, 2, double> N;
2138 * M.random();
2139 * N.random();
2140 * auto S = N * M; // Results in a 3 x 3 Matrix since N.rows() = 3 and M.columns = 3
2141 * \endcode
2142 *
2143 * __RowVector * #Vector / #Vector * #RowVector:__
2144 *
2145 * Defined as standard [dot product](https://en.wikipedia.org/wiki/Dot_product) between vectors as
2146 * long as vectors are of same length
2147 *
2148 * \code {.cpp}
2149 * RowVector<n,MyType> V;
2150 * Vector<n,MyType> W;
2151 * //
2152 * // Fill V and W
2153 * //
2154 * auto S = V*W; // Results in MyType scalar
2155 * \endcode
2156 *
2157 * #Vector * RowVector is same as outer product which is equivalent to a matrix
2158 * multiplication
2159 *
2160 * \code {.cpp}
2161 * auto S = W*V // Result in n x n Matrix of type MyType
2162 * \endcode
2163 *
2164 * __Matrix * Scalar / Scalar * Matrix:__
2165 *
2166 * Multiplication operator between Matrix and Scalar are defined in the usual way (element wise)
2167 *
2168 * \code {.cpp}
2169 * Matrix<n,n,MyType> M;
2170 * M.fill(1);
2171 * auto S = M*2; // Resulting Matrix is uniformly 2
2172 * \endcode
2173 *
2174 * @tparam Mt1 Matrix type for a
2175 * @tparam Mt2 Matrix type for b
2176 * @tparam n Number of rows
2177 * @tparam m Number of columns
2178 * @param a Left Matrix, Vector or Scalar
2179 * @param b Right Matrix, Vector or Scalar
2180 * @return Matrix<n, m, R>
2181 */
2182template <typename Mt1, typename Mt2,
2183 std::enable_if_t<Mt1::is_matrix() && Mt2::is_matrix() && !std::is_same<Mt1, Mt2>::value,
2184 int> = 0,
2185 typename R = hila::type_mul<hila::number_type<Mt1>, hila::number_type<Mt2>>,
2186 int n = Mt1::rows(), int m = Mt2::columns()>
2187inline Matrix<n, m, R> operator*(const Mt1 &a, const Mt2 &b) {
2188
2189 constexpr int p = Mt1::columns();
2190 static_assert(p == Mt2::rows(), "Matrix size: LHS columns != RHS rows");
2191
2192 Matrix<n, m, R> res;
2193
2194 if constexpr (n > 1 && m > 1) {
2195 for (int i = 0; i < n; i++)
2196 for (int j = 0; j < m; j++) {
2197 res.e(i, j) = a.e(i, 0) * b.e(0, j);
2198 for (int k = 1; k < p; k++) {
2199 res.e(i, j) += a.e(i, k) * b.e(k, j);
2200 }
2201 }
2202 } else if constexpr (m == 1) {
2203 // matrix x vector
2204 for (int i = 0; i < n; i++) {
2205 res.e(i) = a.e(i, 0) * b.e(0);
2206 for (int k = 1; k < p; k++) {
2207 res.e(i) += a.e(i, k) * b.e(k);
2208 }
2209 }
2210 } else if constexpr (n == 1) {
2211 // vector x matrix
2212 for (int j = 0; j < m; j++) {
2213 res.e(j) = a.e(0) * b.e(0, j);
2214 for (int k = 1; k < p; k++) {
2215 res.e(j) += a.e(k) * b.e(k, j);
2216 }
2217 }
2218 }
2219
2220 return res;
2221}
2222
2223// and treat separately horiz. vector * vector
2224template <int m, int n, typename T1, typename T2, typename Rt = hila::ntype_op<T1, T2>>
2225Rt operator*(const Matrix<1, m, T1> &A, const Matrix<n, 1, T2> &B) {
2226
2227 static_assert(m == n, "Vector lengths do not match");
2228
2229 Rt res;
2230 res = A.e(0) * B.e(0);
2231 for (int i = 1; i < m; i++) {
2232 res += A.e(i) * B.e(i);
2233 }
2234 return res;
2235}
2236
2237///////////////////////////////////////////////////////////////////////////
2238
2239// matrix * scalar
2240template <typename Mt, typename S,
2241 std::enable_if_t<(Mt::is_matrix() && hila::is_complex_or_arithmetic<S>::value), int> = 0,
2242 typename Rt = hila::mat_scalar_type<Mt, S>>
2243inline Rt operator*(const Mt &mat, const S &rhs) {
2244 Rt res;
2245 for (int i = 0; i < Rt::rows() * Rt::columns(); i++) {
2246 res.c[i] = mat.c[i] * rhs;
2247 }
2248 return res;
2249}
2250
2251// scalar * matrix
2252template <typename Mt, typename S,
2253 std::enable_if_t<(Mt::is_matrix() && hila::is_complex_or_arithmetic<S>::value), int> = 0,
2254 typename Rt = hila::mat_scalar_type<Mt, S>>
2255inline Rt operator*(const S &rhs, const Mt &mat) {
2256 return mat * rhs; // assume commutes
2257}
2258
2259// matrix / scalar
2260
2261/**
2262 * @brief Division operator
2263 *
2264 * Defined for following operations
2265 *
2266 * __ Matrix / Scalar: __
2267 *
2268 * Division operator between Matrix and Scalar are defined in the usual way (element wise)
2269 *
2270 * \code {.cpp}
2271 * Matrix<n,m,MyType> M;
2272 * M.fill(2);
2273 * auto S = M/2; // Resulting matrix is uniformly 1
2274 * \endcode
2275 *
2276 * @tparam Mt Matrix type
2277 * @tparam S Scalar type
2278 * @param mat Matrix to divide scalar with
2279 * @param rhs Scalar to divide with
2280 * @return Rt Resulting Matrix
2281 */
2282template <typename Mt, typename S,
2283 std::enable_if_t<(Mt::is_matrix() && hila::is_complex_or_arithmetic<S>::value), int> = 0,
2284 typename Rt = hila::mat_scalar_type<Mt, S>>
2285inline Rt operator/(const Mt &mat, const S &rhs) {
2286 Rt res;
2287 for (int i = 0; i < Rt::rows() * Rt::columns(); i++) {
2288 res.c[i] = mat.c[i] / rhs;
2289 }
2290 return res;
2291}
2292
2293
2294/**
2295 * @brief Returns trace of product of two matrices
2296 * @details Wrapper around Matrix::mul_trace in the form
2297 * \code {.cpp}
2298 * mul_trace(a,b) // -> a.mul_trace(b)
2299 * \endcode
2300 *
2301 * @tparam Mtype1 Matrix type for a
2302 * @tparam Mtype2 Matrix type for b
2303 * @param a Left Matrix
2304 * @param b Right Matrix
2305 * @return auto Resulting trace
2306 */
2307template <typename Mtype1, typename Mtype2,
2308 std::enable_if_t<Mtype1::is_matrix() && Mtype2::is_matrix(), int> = 0>
2309inline auto mul_trace(const Mtype1 &a, const Mtype2 &b) {
2310
2311 static_assert(Mtype1::columns() == Mtype2::rows() && Mtype1::rows() == Mtype2::columns(),
2312 "mul_trace(a,b): matrix sizes are not compatible");
2313 return a.mul_trace(b);
2314}
2315
2316/**
2317 * @brief compute product of two square matrices and write result to existing matrix
2318 *
2319 * @tparam Mt Matrix type
2320 * @param a Left Matrix
2321 * @param b Right Matrix
2322 * @param c Matrix to which result gets written
2323 * @return void
2324 */
2325// same type square matrices:
2326template <typename Mt, std::enable_if_t<Mt::is_matrix() && Mt::is_square(), int> = 0>
2327inline void mult(const Mt &a, const Mt &b, out_only Mt &res) {
2328 constexpr int n = Mt::rows();
2329 int i, j, k;
2330 for (i = 0; i < n; ++i) {
2331 for (j = 0; j < n; ++j) {
2332 res.e(i, j) = a.e(i, 0) * b.e(0, j);
2333 for (k = 1; k < n; ++k) {
2334 res.e(i, j) += a.e(i, k) * b.e(k, j);
2335 }
2336 }
2337 }
2338}
2339
2340//////////////////////////////////////////////////////////////////////////////////
2341
2342
2343// Stream operator
2344/**
2345 * @brief Stream operator
2346 * @details Naive Stream operator for formatting Matrix for printing. Simply puts elements one after
2347 * the other in row major order
2348 *
2349 * @tparam n
2350 * @tparam m
2351 * @tparam T
2352 * @tparam MT
2353 * @param strm
2354 * @param A
2355 * @return std::ostream&
2356 */
2357template <int n, int m, typename T, typename MT>
2358std::ostream &operator<<(std::ostream &strm, const Matrix_t<n, m, T, MT> &A) {
2359 for (int i = 0; i < n; i++)
2360 for (int j = 0; j < m; j++) {
2361 strm << A.e(i, j);
2362 if (i < n - 1 || j < m - 1)
2363 strm << ' ';
2364 }
2365 return strm;
2366}
2367
2368// Convert to string for "pretty" printing
2369//
2370
2371namespace hila {
2372
2373/**
2374 * @brief Converts Matrix_t object to string
2375 *
2376 * @tparam n Number of rows
2377 * @tparam m Number of columns
2378 * @tparam T Matrix element type
2379 * @tparam MT Matrix type
2380 * @param A Matrix to convert to string
2381 * @param prec Precision of T
2382 * @param separator Separator between elements
2383 * @return std::string
2384 */
2385template <int n, int m, typename T, typename MT>
2386std::string to_string(const Matrix_t<n, m, T, MT> &A, int prec = 8, char separator = ' ') {
2387 std::stringstream strm;
2388 strm.precision(prec);
2389
2390 for (int i = 0; i < n; i++)
2391 for (int j = 0; j < m; j++) {
2392 strm << hila::to_string(A.e(i, j), prec, separator);
2393 if (i < n - 1 || j < m - 1)
2394 strm << separator;
2395 }
2396 return strm.str();
2397}
2398
2399/**
2400 * @brief Formats Matrix_t object in a human readable way
2401 * @details Example 2 x 3 matrix is the following
2402 * \code {.cpp}
2403 * Matrix<2, 3, double> W;
2404 * W.random();
2405 * hila::out0 << hila::prettyprint(W ,4) << '\n';
2406 * // Result:
2407 * // [ 0.8555 0.006359 0.3193 ]
2408 * // [ 0.237 0.8871 0.8545 ]
2409 * \endcode
2410 *
2411 *
2412 * @tparam n Number of rows
2413 * @tparam m Number of columns
2414 * @tparam T Matrix element type
2415 * @tparam MT Matrix type
2416 * @param A Matrix to be formatted
2417 * @param prec Precision of Matrix element
2418 * @return std::string
2419 */
2420template <int n, int m, typename T, typename MT>
2421std::string prettyprint(const Matrix_t<n, m, T, MT> &A, int prec = 8) {
2422 std::stringstream strm;
2423 strm.precision(prec);
2424
2425 if constexpr (n == 1) {
2426 // print a vector, horizontally
2427 strm << '[';
2428 for (int i = 0; i < n * m; i++)
2429 strm << ' ' << prettyprint(A.e(i), prec);
2430 strm << " ]";
2431 // if (n > 1)
2432 // strm << "^T";
2433 } else {
2434 // now a matrix - split the matrix on lines.
2435 // do it so that columns are equal width...
2436 std::vector<std::string> lines, columns;
2437 lines.resize(n);
2438 columns.resize(n);
2439
2440 for (int i = 0; i < n; i++)
2441 lines[i] = "[ ";
2442 for (int j = 0; j < m; j++) {
2443 int size = 0;
2444 for (int i = 0; i < n; i++) {
2445 std::stringstream item;
2446 item << prettyprint(A.e(i, j), prec);
2447 columns[i] = item.str();
2448 if (columns[i].size() > size)
2449 size = columns[i].size();
2450 }
2451 for (int i = 0; i < n; i++) {
2452 lines[i].append(size - columns[i].size(), ' ');
2453 lines[i].append(columns[i]);
2454 lines[i].append(1, ' ');
2455 }
2456 }
2457 for (int i = 0; i < n - 1; i++) {
2458 strm << lines[i] << "]\n";
2459 }
2460 strm << lines[n - 1] << "]";
2461 }
2462 return strm.str();
2463}
2464
2465} // namespace hila
2466
2467/**
2468 * @brief Returns square norm of Matrix
2469 * @details Wrapper around Matrix::squarenorm - sum of squared elements
2470 *
2471 * Can be called as:
2472 * \code {.cpp}
2473 * Matrix<n,m,MyType> M;
2474 * auto a = squarenorm(M);
2475 * \endcode
2476 *
2477 * @tparam Mt Matrix type
2478 * @param rhs Matrix to compute square norm of
2479 * @return auto
2480 */
2481template <typename Mt, std::enable_if_t<Mt::is_matrix(), int> = 0>
2482inline auto squarenorm(const Mt &rhs) {
2483 return rhs.squarenorm();
2484}
2485
2486// Vector norm - sqrt of squarenorm()
2487
2488/**
2489 * @brief Returns vector norm of Matrix
2490 * @details Wrapper around Matrix::norm - sqrt of sum of squared elements
2491 * @tparam Mt Matrix type
2492 * @param rhs Matrix to compute norm of
2493 * @return auto
2494 */
2495template <typename Mt, std::enable_if_t<Mt::is_matrix(), int> = 0>
2496inline auto norm(const Mt &rhs) {
2497 return rhs.norm();
2498}
2499
2500
2501// Cast operators to different number type
2502// cast_to<double>(a);
2503
2504template <typename Ntype, typename T, int n, int m,
2505 std::enable_if_t<hila::is_arithmetic<T>::value, int> = 0>
2508 for (int i = 0; i < n * m; i++)
2509 res.c[i] = mat.c[i];
2510 return res;
2511}
2512
2513template <typename Ntype, typename T, int n, int m,
2514 std::enable_if_t<hila::is_complex<T>::value, int> = 0>
2517 for (int i = 0; i < n * m; i++)
2518 res.c[i] = cast_to<Ntype>(mat.c[i]);
2519 return res;
2520}
2521
2522// Calculate exp of a square matrix
2523// Go to order ORDER in the exponential of the matrices
2524// matrices, since unitarity seems important.
2525// Evaluation is done as:
2526// exp(H) = 1 + H + H^2/2! + H^3/3! ..-
2527// = 1 + H*( 1 + (H/2)*( 1 + (H/3)*( ... )))
2528// Do it backwards in order to reduce accumulation of errors
2529
2530/**
2531 * @brief Calculate exp of a square matrix
2532 * @details Computes up to certain order given as argument
2533 *
2534 * __Evaluation is done as__:
2535 * \f{align}{ \exp(H) &= 1 + H + \frac{H^2}{2!} + \frac{H^2}{2!} + \frac{H^3}{3!} \\
2536 * &= 1 + H\cdot(1 + (\frac{H}{2})\cdot (1 + (\frac{H}{3})\cdot(...))) \f}
2537 * Done backwards in order to reduce accumulation of errors
2538
2539 * @tparam n Number of rowsMa
2540 * @tparam T Matrix element type
2541 * @tparam MT Matrix type
2542 * @param mat Matrix to compute exponential for
2543 * @param order order to compute exponential to
2544 * @return Matrix_t<n, m, T, MT>
2545 */
2546template <int n, int m, typename T, typename MT>
2547inline Matrix_t<n, m, T, MT> exp(const Matrix_t<n, m, T, MT> &mat, const int order = 20) {
2548 static_assert(n == m, "exp() only for square matrices");
2549
2550 hila::arithmetic_type<T> one = 1.0;
2551 Matrix_t<n, m, T, MT> r = mat;
2552
2553 // doing the loop step-by-step should reduce temporaries
2554 for (int k = order; k > 1; k--) {
2555 r *= (one / k);
2556 r += one;
2557 r *= mat;
2558 }
2559 r += one;
2560
2561 return r;
2562}
2563
2564// Calculate exp of a square matrix
2565// using iterative Cayley-Hamilton described in arXiv:2404.07704
2566/**
2567 * @brief Calculate exp of a square matrix
2568 * @details Computation is done using iterative Cayley-Hamilton (cf. from arXiv:2404.07704)
2569
2570 * @tparam n Number of rowsMa
2571 * @tparam T Matrix element type
2572 * @tparam MT Matrix type
2573 * @param mat Matrix to compute exponential for
2574 * @param omat Matrix to which exponential of mat gets stored (optional)
2575 * @param pl array of n+1 temporary nxn Matrices (optional)
2576 * @return void (if omat is provided) or Matrix_t<n,m,T,MT>
2577 */
2578template <int n, int m, typename T, typename MT>
2579inline void chexp(const Matrix_t<n, m, T, MT> &mat, out_only Matrix_t<n, m, T, MT> &omat,
2580 Matrix_t<n, m, T, MT>(out_only &pl)[n + 1]) {
2581 static_assert(n == m, "chexp() only for square matrices");
2582 // compute the first n matrix powers of mat and the corresponding traces :
2583 // the i-th matrix power of mat[][] is stored in pl[i][][]
2584 T trpl[n + 1]; // the trace of pl[i][][] is stored in trpl[i]
2585 trpl[0] = n;
2586 pl[1] = mat;
2587 trpl[1] = trace(pl[1]);
2588 int i, j, k;
2589 for (i = 2; i <= n; ++i) {
2590 j = i / 2;
2591 k = i % 2;
2592 mult(pl[j], pl[j + k], pl[i]);
2593 trpl[i] = trace(pl[i]);
2594 }
2595
2596 // compute the characteristic polynomial coefficients crpl[] from the traced powers trpl[] :
2597 T crpl[n + 1];
2598 crpl[n] = 1;
2599 for (j = 1; j <= n; ++j) {
2600 crpl[n - j] = 0;
2601 for (i = 1; i <= j; ++i) {
2602 crpl[n - j] -= crpl[n - (j - i)] * trpl[i];
2603 }
2604 crpl[n - j] /= j;
2605 }
2606
2607
2608 int mmax = 15 * n; // maximum number of Cayley-Hamilton iterations if no convergence is reached
2609 T al[n], pal[n]; // temp. Cayley-Hamilton coefficents
2610 hila::arithmetic_type<T>
2611 wpf = 1.0,
2612 twpf = 1.0, ttwpf; // initial values for power series coefficnet and its running sum
2613
2614 // set initial values for the n entries in al[] and pal[] :
2615 for (i = 0; i < n; ++i) {
2616 pal[i] = 0;
2617 al[i] = wpf;
2618 wpf /= (i + 1); // compute (i+1)-th power series coefficent from the i-th coefficient
2619 twpf += wpf;
2620 }
2621 pal[n - 1] = 1.0;
2622
2623 // next we iteratively add higher order power series terms to al[] till al[] stops changing
2624 // more precisely: the iteration will terminate as soon as twpf stops changing. Here twpf
2625 // is the sum \sum_{i=0}^{j} s_i/i!, with s_i referring to the magnitude the vector pal[]
2626 // would have at iteration i, if no renormalization were used.
2627 T ch, cho; // temporary variables for iteration
2628 hila::arithmetic_type<T> s, rs = 1.0; // temp variables used for renormalization of pal[]
2629 ttwpf = twpf;
2630 for (j = n; j < mmax; ++j) {
2631 pal[n - 1] *= rs;
2632 ch = -pal[n - 1] * crpl[0];
2633 cho = pal[0] * rs;
2634 pal[0] = ch;
2635 s = squarenorm(ch);
2636 al[0] += wpf * ch;
2637 for (i = 1; i < n; ++i) {
2638 ch = cho - pal[n - 1] * crpl[i];
2639 cho = pal[i] * rs;
2640 pal[i] = ch;
2641 s += squarenorm(ch);
2642 al[i] += wpf * ch;
2643 }
2644
2645 if (s > 1.0) {
2646 // if s is bigger than 1, normalize pal[] by a factor rs=1.0/s in next itaration,
2647 // and multiply wpf by s to compensate
2648 s = std::sqrt(s);
2649 wpf *= s / (j + 1);
2650 rs = 1.0 / s;
2651 } else {
2652 wpf /= (j + 1);
2653 rs = 1.0;
2654 }
2655 twpf += wpf;
2656 if (ttwpf == twpf) {
2657 // terminate iteration when numeric value of twpf stops changing
2658 break;
2659 }
2660 ttwpf = twpf;
2661 }
2662 // if(hila::myrank()==0) {
2663 // std::cout<<"chexp niter: "<<j<<" ("<<j-n<<")"<<std::endl;
2664 // }
2665
2666 // form output matrix:
2667 for (i = 0; i < n; ++i) {
2668 for (j = 0; j < n; ++j) {
2669 if (i == j) {
2670 omat.e(i, j) = al[0];
2671 } else {
2672 omat.e(i, j) = 0;
2673 }
2674 for (k = 1; k < n; ++k) {
2675 omat.e(i, j) += al[k] * pl[k].e(i, j);
2676 }
2677 }
2678 }
2679}
2680
2681// overload wrapper for chexp where omat is not provided
2682template <int n, int m, typename T, typename MT>
2684 Matrix_t<n, m, T, MT>(out_only &pl)[n + 1]) {
2685 static_assert(n == m, "chexp() only for square matrices");
2686 chexp(mat, pl[0], pl);
2687 return pl[0];
2688}
2689
2690// overload wrapper for chexp which creates the temporary matrix array pl[n+1] internally
2691template <int n, int m, typename T, typename MT>
2692inline void chexp(const Matrix_t<n, m, T, MT> &mat, out_only Matrix_t<n, m, T, MT> &omat) {
2693 static_assert(n == m, "chexp() only for square matrices");
2694 Matrix_t<n, m, T, MT> pl[n + 1];
2695 chexp(mat, omat, pl);
2696}
2697
2698// overload wrapper for chexp where omat is not provided
2699// and which creates the temporary matrix array pl[n+1] internally
2700template <int n, int m, typename T, typename MT>
2702 static_assert(n == m, "chexp() only for square matrices");
2703 Matrix_t<n, m, T, MT> pl[n + 1];
2704 chexp(mat, pl[0], pl);
2705 return pl[0];
2706}
2707
2708
2709// Calculate exp of a square matrix
2710// using iterative Cayley-Hamilton described in arXiv:2404.07704
2711/**
2712 * @brief Calculate exp of a square matrix
2713 * @details Computation is done using iterative Cayley-Hamilton (cf. from arXiv:2404.07704) with
2714 minimal temporary storage
2715
2716 * @tparam n Number of rowsMa
2717 * @tparam T Matrix element type
2718 * @tparam MT Matrix type
2719 * @param mat Matrix to compute exponential for
2720 * @return Matrix_t<n, m, T, MT>
2721 */
2722template <int n, int m, typename T, typename MT>
2724 static_assert(n == m, "chsexp() only for square matrices");
2725
2726 // compute the characteristic polynomial coefficients crpl[] with the Faddeev-LeVerrier
2727 // algorithm :
2728 int i, j, k;
2730 T crpl[n + 1];
2731 crpl[n] = 1.0;
2732 int ip = 0;
2733 tB[ip] = 1.;
2734 T tc = trace(mat);
2735 crpl[n - 1] = tc;
2736 tB[1 - ip] = mat;
2737 for (k = 2; k <= n; ++k) {
2738 tB[1 - ip] -= tc;
2739 mult(mat, tB[1 - ip], tB[ip]);
2740 tc = trace(tB[ip]) / k;
2741 crpl[n - k] = tc;
2742 ip = 1 - ip;
2743 }
2744
2745
2746 int mmax = 15 * n; // maximum number of Cayley-Hamilton iterations if no convergence is reached
2747 T al[n], pal[n]; // temp. Cayley-Hamilton coefficents
2748 hila::arithmetic_type<T> wpf = 1.0, twpf = 1.0, ttwpf; // leading coefficient of power series
2749 // set initial values for the n entries in al[] and pal[] :
2750 for (i = 0; i < n; ++i) {
2751 pal[i] = 0;
2752 al[i] = wpf;
2753 wpf /= (i + 1); // compute (i+1)-th power series coefficent from the i-th coefficient
2754 twpf += wpf;
2755 }
2756 pal[n - 1] = 1.0;
2757
2758 // next we iteratively add higher order power series terms to al[] till al[] stops changing
2759 // more precisely: the iteration will terminate as soon as twpf stops changing. Here twpf
2760 // is the sum \sum_{i=0}^{j} s_i/i!, with s_i referring to the magnitude the vector pal[]
2761 // would have at iteration i, if no renormalization were used.
2762 T ch, cho; // temporary variables for iteration
2763 hila::arithmetic_type<T> s, rs = 1.0; // temp variables used for renormalization of pal[]
2764 ttwpf = twpf;
2765 for (j = n; j < mmax; ++j) {
2766 pal[n - 1] *= rs;
2767 ch = -pal[n - 1] * crpl[0];
2768 cho = pal[0] * rs;
2769 pal[0] = ch;
2770 s = squarenorm(ch);
2771 al[0] += wpf * ch;
2772 for (i = 1; i < n; ++i) {
2773 ch = cho - pal[n - 1] * crpl[i];
2774 cho = pal[i] * rs;
2775 pal[i] = ch;
2776 s += squarenorm(ch);
2777 al[i] += wpf * ch;
2778 }
2779 if (s > 1.0) {
2780 // if s is bigger than 1, normalize pal[] by a factor rs=1.0/s in next itaration,
2781 // and multiply wpf by s to compensate
2782 s = std::sqrt(s);
2783 wpf *= s / (j + 1);
2784 rs = 1.0 / s;
2785 } else {
2786 wpf /= (j + 1);
2787 rs = 1.0;
2788 }
2789 twpf += wpf;
2790 if (ttwpf == twpf) {
2791 // terminate iteration
2792 break;
2793 }
2794 ttwpf = twpf;
2795 }
2796 // if(hila::myrank()==0) {
2797 // std::cout<<"chsexp niter: "<<j<<" ("<<j-n<<")"<<std::endl;
2798 // }
2799
2800 // form output matrix:
2801 ip = 0;
2802 tB[ip] = mat;
2803 tB[ip] *= al[n - 1];
2804 tB[ip] += al[n - 2];
2805 for (i = 2; i < n; ++i) {
2806 mult(tB[ip], mat, tB[1 - ip]);
2807 tB[1 - ip] += al[n - i - 1];
2808 ip = 1 - ip;
2809 }
2810
2811 return tB[ip];
2812}
2813
2814
2815#include "datatypes/array.h"
2816
2817#include "datatypes/diagonal_matrix.h"
2818
2819#include "datatypes/matrix_linalg.h"
2820
2821// #include "datatypes/dagger.h"
2822
2823#endif
Definition of Array class.
Array< n, m, T > sqrt(Array< n, m, T > a)
Square root.
Definition array.h:1039
Array< n, m, T > operator*(Array< n, m, T > a, const Array< n, m, T > &b)
Multiplication operator.
Definition array.h:861
Array< n, m, Ntype > cast_to(const Array< n, m, T > &mat)
Array casting operation.
Definition array.h:1289
Array type
Definition array.h:43
Complex definition.
Definition cmplx.h:56
Define type DiagonalMatrix<n,T>
T e(const int i) const
Element access - e(i) gives diagonal element i.
The main matrix type template Matrix_t. This is a base class type for "useful" types which are deriv...
Definition matrix.h:106
void set_column(int c, const Vector< n, S > &v)
get column of a matrix
Definition matrix.h:395
void mult_by_2x2_left(int p, int q, const Matrix_t< 2, 2, R, Mt > &M)
Multiply -matrix from the left by matrix defined by sub matrix.
Definition matrix.h:1473
T trace() const
Computes Trace for Matrix.
Definition matrix.h:1079
Vector< n, T > column(int c) const
Returns column vector as value at index c.
Definition matrix.h:368
hila::arithmetic_type< T > norm() const
Calculate vector norm - sqrt of squarenorm.
Definition matrix.h:1133
void set_row(int r, const RowVector< m, S > &v)
Set row of Matrix with RowVector if types are assignable.
Definition matrix.h:357
auto & operator=(const Matrix_t< n, m, S, MT > &rhs) &
Assignment operator = to assign values to matrix.
Definition matrix.h:594
static constexpr int rows()
Define constant methods rows(), columns() - may be useful in template code.
Definition matrix.h:226
const auto & fill(const S rhs)
Matrix fill.
Definition matrix.h:937
static constexpr int columns()
Returns column length.
Definition matrix.h:234
Mtype & random()
dot with matrix - matrix
Definition matrix.h:1296
void mult_by_2x2_right(int p, int q, const Matrix_t< 2, 2, R, Mt > &M)
Multiply -matrix from the right by matrix defined by sub matrix.
Definition matrix.h:1502
Mtype sort(Vector< N, int > &permutation, hila::sort order=hila::sort::ascending) const
Sort method for Vector.
Definition matrix.h:1416
Matrix< n, m, hila::arithmetic_type< T > > imag() const
Returns imaginary part of Matrix or Vector.
Definition matrix.h:1064
void set_diagonal(const Vector< n, S > &v)
Set diagonal of square matrix to Vector which is passed to the method.
Definition matrix.h:430
auto & operator-=(const Matrix_t< n, m, S, MT > &rhs) &
Subtract assign operator with matrix or scalar.
Definition matrix.h:726
const RowVector< n, T > & transpose() const
Transpose of vector.
Definition matrix.h:971
const auto & operator+() const
Unary + operator.
Definition matrix.h:507
static constexpr int size()
Returns size of Vector or square Matrix.
Definition matrix.h:248
Mtype permute(const Vector< N, int > &permutation) const
Permute elements of vector.
Definition matrix.h:1377
Matrix_t()=default
Define default constructors to ensure std::is_trivial.
Mtype conj() const
Returns complex conjugate of Matrix.
Definition matrix.h:986
Mtype permute_rows(const Vector< n, int > &permutation) const
permute rows of Matrix
Definition matrix.h:1361
static constexpr bool is_vector()
Returns true if Matrix is a vector.
Definition matrix.h:132
T det() const
determinant function - if matrix size is < 5, use Laplace, otherwise LU
auto abs() const
Returns absolute value of Matrix.
Definition matrix.h:1031
Matrix< n, p, R > outer_product(const Matrix< p, q, S > &rhs) const
Outer product.
Definition matrix.h:1263
T max() const
Find max or min value - only for arithmetic types.
Definition matrix.h:1148
Rtype dagger() const
Hermitian conjugate of matrix.
Definition matrix.h:1002
Matrix< n, m, hila::arithmetic_type< T > > real() const
Returns real part of Matrix or Vector.
Definition matrix.h:1051
bool operator!=(const Matrix< n, m, S > &rhs) const
Boolean operator != to check if matrices are exactly different.
Definition matrix.h:544
int eigen_hermitean(DiagonalMatrix< n, Et > &eigenvalues, Matrix_t< n, n, Mt, MT > &eigenvectors, enum hila::sort sorted=hila::sort::unsorted) const
Calculate eigenvalues and -vectors of an hermitean (or symmetric) matrix.
int svd_pivot(Matrix_t< n, n, Mt, MT > &_U, DiagonalMatrix< n, Et > &_D, Matrix_t< n, n, Mt, MT > &_V, enum hila::sort sorted=hila::sort::unsorted) const
Singular value decomposition: divide matrix A = U S V*, where U,V unitary and S diagonal matrix of re...
bool operator==(const Matrix< n2, m2, S > &rhs) const
Boolean operator == to determine if two matrices are exactly the same.
Definition matrix.h:522
T c[n *m]
The data as a one dimensional array.
Definition matrix.h:110
auto & operator+=(const Matrix_t< n, m, S, MT > &rhs) &
Add assign operator with matrix or scalar.
Definition matrix.h:690
Mtype & gaussian_random(base_type width=1.0)
Fills Matrix with gaussian random elements.
Definition matrix.h:1314
T e(const int i, const int j) const
Standard array indexing operation for matrices and vectors.
Definition matrix.h:286
T det_laplace() const
Determinant of matrix, using Laplace method.
RowVector< m, T > row(int r) const
Return row of a matrix as a RowVector.
Definition matrix.h:342
const DiagonalMatrix< n, T > & asDiagonalMatrix() const
Cast Vector to DiagonalMatrix.
Definition matrix.h:458
Rtype transpose() const
Transpose of matrix.
Definition matrix.h:954
auto & operator*=(const DiagonalMatrix< m, S > &rhs) &
mult and divide assign a diagonal - cols must match diagonal matrix rows
Definition matrix.h:900
hila::type_mul< T, S > mul_trace(const Matrix_t< p, q, S, MT > &rm) const
Multiply with given matrix and compute trace of result.
Definition matrix.h:1100
T operator[](const int i) const
Indexing operation [] defined only for vectors.
Definition matrix.h:327
int svd(Matrix_t< n, n, Mt, MT > &_U, DiagonalMatrix< n, Et > &_D, Matrix_t< n, n, Mt, MT > &_V, enum hila::sort sorted=hila::sort::unsorted) const
Singular value decomposition: divide matrix A = U S V*, where U,V unitary and S diagonal matrix of re...
auto & operator+=(const DiagonalMatrix< n, S > &rhs) &
add and sub assign a DiagonalMatrix
Definition matrix.h:877
Mtype permute_columns(const Vector< m, int > &permutation) const
permute columns of Matrix
Definition matrix.h:1347
R dot(const Matrix< p, q, S > &rhs) const
Dot product.
Definition matrix.h:1220
T det_lu() const
following calculate the determinant of a square matrix det() is the generic interface,...
const Array< n, m, T > & asArray() const
Cast Matrix to Array.
Definition matrix.h:443
Mtype operator-() const
Unary - operator.
Definition matrix.h:493
auto & operator/=(const S rhs) &
Divide assign scalar.
Definition matrix.h:863
DiagonalMatrix< n, T > diagonal()
Return diagonal of square matrix.
Definition matrix.h:406
std::string str(int prec=8, char separator=' ') const
Definition matrix.h:1587
Rtype adjoint() const
Adjoint of matrix.
Definition matrix.h:1019
auto & operator*=(const Matrix_t< m, p, S, MT > &rhs) &
Multiply assign scalar or matrix.
Definition matrix.h:796
hila::arithmetic_type< T > squarenorm() const
Calculate square norm - sum of squared elements.
Definition matrix.h:1117
static constexpr bool is_square()
Returns true if matrix is a square matrix.
Definition matrix.h:142
Matrix class which defines matrix operations.
Definition matrix.h:1620
hila::arithmetic_type< T > base_type
std incantation for field types
Definition matrix.h:1624
Definition of Complex types.
T arg(const Complex< T > &a)
Return argument of Complex number.
Definition cmplx.h:1199
This file defines all includes for HILA.
Matrix_t< n, m, T, MT > chsexp(const Matrix_t< n, m, T, MT > &mat)
Calculate exp of a square matrix.
Definition matrix.h:2723
void chexp(const Matrix_t< n, m, T, MT > &mat, Matrix_t< n, m, T, MT > &omat, Matrix_t< n, m, T, MT >(&pl)[n+1])
Calculate exp of a square matrix.
Definition matrix.h:2579
auto dagger(const Mtype &arg)
Return dagger of Matrix.
Definition matrix.h:1813
Rtype operator+(const Mtype1 &a, const Mtype2 &b)
Addition operator.
Definition matrix.h:1946
auto abs(const Mtype &arg)
Return absolute value Matrix or Vector.
Definition matrix.h:1833
Matrix_t< n, m, T, MT > exp(const Matrix_t< n, m, T, MT > &mat, const int order=20)
Calculate exp of a square matrix.
Definition matrix.h:2547
auto norm(const Mt &rhs)
Returns vector norm of Matrix.
Definition matrix.h:2496
auto conj(const Mtype &arg)
Return conjugate Matrix or Vector.
Definition matrix.h:1777
auto squarenorm(const Mt &rhs)
Returns square norm of Matrix.
Definition matrix.h:2482
auto transpose(const Mtype &arg)
Return transposed Matrix or Vector.
Definition matrix.h:1757
Rtype operator-(const Mtype1 &a, const Mtype2 &b)
Subtraction operator.
Definition matrix.h:2022
auto real(const Mtype &arg)
Return real of Matrix or Vector.
Definition matrix.h:1874
auto imag(const Mtype &arg)
Return imaginary of Matrix or Vector.
Definition matrix.h:1894
std::ostream & operator<<(std::ostream &strm, const Matrix_t< n, m, T, MT > &A)
Stream operator.
Definition matrix.h:2358
auto trace(const Mtype &arg)
Return trace of square Matrix.
Definition matrix.h:1854
auto mul_trace(const Mtype1 &a, const Mtype2 &b)
Returns trace of product of two matrices.
Definition matrix.h:2309
void mult(const Mt &a, const Mt &b, Mt &res)
compute product of two square matrices and write result to existing matrix
Definition matrix.h:2327
auto adjoint(const Mtype &arg)
Return adjoint Matrix.
Definition matrix.h:1797
Implement hila::swap for gauge fields.
Definition array.h:981
double random()
Real valued uniform random number generator.
Definition hila_gpu.cpp:118
double gaussrand()
Gaussian random generation routine.
Definition random.cpp:216
double gaussrand2(double &out2)
hila::gaussrand2 returns 2 Gaussian distributed random numbers with variance .
Definition random.cpp:181
decltype(std::declval< A >()+std::declval< B >()) type_plus
Definition type_tools.h:103
std::string to_string(const Array< n, m, T > &A, int prec=8, char separator=' ')
Converts Array object to string.
Definition array.h:995
type to store the return value of eigen_hermitean(): {E, U} where E is nxn DiagonalMatrix containing ...
Definition matrix.h:73
hila::is_complex_or_arithmetic<T>::value
Definition cmplx.h:665
type to store the return combo of svd: {U, D, V} where U and V are nxn unitary / orthogonal,...
Definition matrix.h:61