HILA
Loading...
Searching...
No Matches
matrix.h File Reference

Definition of Matrix types. More...

#include <type_traits>
#include <sstream>
#include "plumbing/defs.h"
#include "datatypes/cmplx.h"
#include "datatypes/array.h"
#include "datatypes/diagonal_matrix.h"
#include "datatypes/matrix_linalg.h"
Include dependency graph for matrix.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

struct  svd_result< M >
 type to store the return combo of svd: {U, D, V} where U and V are nxn unitary / orthogonal, and D is real diagonal singular value matrices. More...
 
struct  eigen_result< M >
 type to store the return value of eigen_hermitean(): {E, U} where E is nxn DiagonalMatrix containing eigenvalues and U nxn unitary matrix, with eigenvector columns More...
 
class  Matrix_t< n, m, T, Mtype >
 The main \( n \times m \) matrix type template Matrix_t. This is a base class type for "useful" types which are derived from this. More...
 
class  Matrix< n, m, T >
 \( n \times m \) Matrix class which defines matrix operations. More...
 

Namespaces

namespace  hila
 Implement hila::swap for gauge fields.
 

Typedefs

template<int n, typename T >
using Vector = Matrix< n, 1, T >
 Vector is defined as 1-column Matrix.
 
template<int n, typename T >
using RowVector = Matrix< 1, n, T >
 RowVector is a 1-row Matrix.
 
template<int n, typename T >
using SquareMatrix = Matrix< n, n, T >
 Square matrix is defined as alias with special case of Matrix<n,n,T>
 

Functions

template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto transpose (const Mtype &arg)
 Return transposed Matrix or Vector.
 
template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto conj (const Mtype &arg)
 Return conjugate Matrix or Vector.
 
template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto adjoint (const Mtype &arg)
 Return adjoint Matrix.
 
template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto dagger (const Mtype &arg)
 Return dagger of Matrix.
 
template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto abs (const Mtype &arg)
 Return absolute value Matrix or Vector.
 
template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto trace (const Mtype &arg)
 Return trace of square Matrix.
 
template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto real (const Mtype &arg)
 Return real of Matrix or Vector.
 
template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto imag (const Mtype &arg)
 Return imaginary of Matrix or Vector.
 
template<typename Mtype1 , typename Mtype2 , std::enable_if_t< Mtype1::is_matrix() &&Mtype2::is_matrix(), int > = 0, typename Rtype = hila::mat_x_mat_type<Mtype1, Mtype2>, std::enable_if_t<!std::is_same< Mtype1, Rtype >::value, int > = 0>
Rtype operator+ (const Mtype1 &a, const Mtype2 &b)
 Addition operator.
 
template<typename Mtype1 , typename Mtype2 , std::enable_if_t< Mtype1::is_matrix() &&Mtype2::is_matrix(), int > = 0, typename Rtype = hila::mat_x_mat_type<Mtype1, Mtype2>>
Rtype operator- (const Mtype1 &a, const Mtype2 &b)
 Subtraction operator.
 
template<typename Mt1 , typename Mt2 , std::enable_if_t< Mt1::is_matrix() &&Mt2::is_matrix() &&!std::is_same< Mt1, Mt2 >::value, int > = 0, typename R = hila::type_mul<hila::number_type<Mt1>, hila::number_type<Mt2>>, int n = Mt1::rows(), int m = Mt2::columns()>
Matrix< n, m, R > operator* (const Mt1 &a, const Mt2 &b)
 Multiplication operator.
 
template<typename Mt , typename S , std::enable_if_t<(Mt::is_matrix() &&hila::is_complex_or_arithmetic< S >::value), int > = 0, typename Rt = hila::mat_scalar_type<Mt, S>>
Rt operator/ (const Mt &mat, const S &rhs)
 Division operator.
 
template<typename Mtype1 , typename Mtype2 , std::enable_if_t< Mtype1::is_matrix() &&Mtype2::is_matrix(), int > = 0>
auto mul_trace (const Mtype1 &a, const Mtype2 &b)
 Returns trace of product of two matrices.
 
template<typename Mt , std::enable_if_t< Mt::is_matrix() &&Mt::is_square(), int > = 0>
void mult (const Mt &a, const Mt &b, Mt &res)
 compute product of two square matrices and write result to existing matrix
 
template<int n, int m, typename T , typename MT >
std::ostream & operator<< (std::ostream &strm, const Matrix_t< n, m, T, MT > &A)
 Stream operator.
 
template<int n, int m, typename T , typename MT >
std::string hila::to_string (const Matrix_t< n, m, T, MT > &A, int prec=8, char separator=' ')
 Converts Matrix_t object to string.
 
template<int n, int m, typename T , typename MT >
std::string hila::prettyprint (const Matrix_t< n, m, T, MT > &A, int prec=8)
 Formats Matrix_t object in a human readable way.
 
template<typename Mt , std::enable_if_t< Mt::is_matrix(), int > = 0>
auto squarenorm (const Mt &rhs)
 Returns square norm of Matrix.
 
template<typename Mt , std::enable_if_t< Mt::is_matrix(), int > = 0>
auto norm (const Mt &rhs)
 Returns vector norm of Matrix.
 
template<int n, int m, typename T , typename MT >
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.
 
template<int n, int m, typename T , typename MT >
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.
 
template<int n, int m, typename T , typename MT >
Matrix_t< n, m, T, MT > chsexp (const Matrix_t< n, m, T, MT > &mat)
 Calculate exp of a square matrix.
 

Detailed Description

Definition of Matrix types.

This file contains base matrix type Matrix_t which defines all general matrix type operations Matrix types are Matrix, Vector, RowVector, SquareMatrix of which Matrix is defined as a class and the rest are special cases of the Matrix class.

Definition in file matrix.h.

Function Documentation

◆ abs()

template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto abs ( const Mtype arg)
inline

Return absolute value Matrix or Vector.

Wrapper around Matrix::abs

Can be called as:

M_abs = abs(M);
Matrix class which defines matrix operations.
Definition matrix.h:1620
T abs(const Complex< T > &a)
Return absolute value of Complex number.
Definition cmplx.h:1187
Template Parameters
MtypeMatrix type
Parameters
argObject to compute absolute value of
Returns
auto Mtype absolute value Matrix or Vector

Definition at line 1833 of file matrix.h.

◆ adjoint()

template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto adjoint ( const Mtype arg)
inline

Return adjoint Matrix.

Wrapper around Matrix::adjoint

Can be called as:

Matrix<n,m,MyType> M,M_adjoint;
M_conjugate = adjoint(M);
Template Parameters
MtypeMatrix type
Parameters
argObject to compute adjoint of
Returns
auto Mtype adjoint matrix

Definition at line 1797 of file matrix.h.

◆ chexp()

template<int n, int m, typename T , typename MT >
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] 
)
inline

Calculate exp of a square matrix.

Computation is done using iterative Cayley-Hamilton (cf. from arXiv:2404.07704)

Template Parameters
nNumber of rowsMa
TMatrix element type
MTMatrix type
Parameters
matMatrix to compute exponential for
omatMatrix to which exponential of mat gets stored (optional)
plarray of n+1 temporary nxn Matrices (optional)
Returns
void (if omat is provided) or Matrix_t<n,m,T,MT>

Definition at line 2579 of file matrix.h.

◆ chsexp()

template<int n, int m, typename T , typename MT >
Matrix_t< n, m, T, MT > chsexp ( const Matrix_t< n, m, T, MT > &  mat)
inline

Calculate exp of a square matrix.

Computation is done using iterative Cayley-Hamilton (cf. from arXiv:2404.07704) with minimal temporary storage

Template Parameters
nNumber of rowsMa
TMatrix element type
MTMatrix type
Parameters
matMatrix to compute exponential for
Returns
Matrix_t<n, m, T, MT>

Definition at line 2723 of file matrix.h.

◆ conj()

template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto conj ( const Mtype arg)
inline

Return conjugate Matrix or Vector.

Wrapper around Matrix::conj

Can be called as:

Matrix<n,m,MyType> M,M_conjugate;
M_conjugate = conj(M);
Array< n, m, T > conj(const Array< n, m, T > &arg)
Return conjugate Array.
Definition array.h:674
Template Parameters
MtypeMatrix type
Parameters
argObject to be conjugated
Returns
auto Mtype conjugated matrix

Definition at line 1777 of file matrix.h.

◆ dagger()

template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto dagger ( const Mtype arg)
inline

Return dagger of Matrix.

Wrapper around Matrix::adjoint

Same as adjoint

Template Parameters
MtypeMatrix type
Parameters
argObject to compute dagger of
Returns
auto Mtype dagger matrix

Definition at line 1813 of file matrix.h.

◆ exp()

template<int n, int m, typename T , typename MT >
Matrix_t< n, m, T, MT > exp ( const Matrix_t< n, m, T, MT > &  mat,
const int  order = 20 
)
inline

Calculate exp of a square matrix.

Computes up to certain order given as argument

Evaluation is done as:

\begin{align} \exp(H) &= 1 + H + \frac{H^2}{2!} + \frac{H^2}{2!} + \frac{H^3}{3!} \\ &= 1 + H\cdot(1 + (\frac{H}{2})\cdot (1 + (\frac{H}{3})\cdot(...))) \end{align}

Done backwards in order to reduce accumulation of errors

Template Parameters
nNumber of rowsMa
TMatrix element type
MTMatrix type
Parameters
matMatrix to compute exponential for
orderorder to compute exponential to
Returns
Matrix_t<n, m, T, MT>

Definition at line 2547 of file matrix.h.

◆ imag()

template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto imag ( const Mtype arg)
inline

Return imaginary of Matrix or Vector.

Wrapper around Matrix::imag

Can be called as:

M_imag = imag(M);
Array< n, m, hila::arithmetic_type< T > > imag(const Array< n, m, T > &arg)
Return imaginary part of Array.
Definition array.h:702
Template Parameters
MtypeMatrix type
Parameters
argObject to compute imag part of
Returns
auto Mtype imag part of arg

Definition at line 1894 of file matrix.h.

◆ mul_trace()

template<typename Mtype1 , typename Mtype2 , std::enable_if_t< Mtype1::is_matrix() &&Mtype2::is_matrix(), int > = 0>
auto mul_trace ( const Mtype1 &  a,
const Mtype2 &  b 
)
inline

Returns trace of product of two matrices.

Wrapper around Matrix::mul_trace in the form

mul_trace(a,b) // -> a.mul_trace(b)
auto mul_trace(const Mtype1 &a, const Mtype2 &b)
Returns trace of product of two matrices.
Definition matrix.h:2309
Template Parameters
Mtype1Matrix type for a
Mtype2Matrix type for b
Parameters
aLeft Matrix
bRight Matrix
Returns
auto Resulting trace

Definition at line 2309 of file matrix.h.

◆ mult()

template<typename Mt , std::enable_if_t< Mt::is_matrix() &&Mt::is_square(), int > = 0>
void mult ( const Mt &  a,
const Mt &  b,
Mt &  res 
)
inline

compute product of two square matrices and write result to existing matrix

Template Parameters
MtMatrix type
Parameters
aLeft Matrix
bRight Matrix
cMatrix to which result gets written
Returns
void

Definition at line 2327 of file matrix.h.

◆ norm()

template<typename Mt , std::enable_if_t< Mt::is_matrix(), int > = 0>
auto norm ( const Mt &  rhs)
inline

Returns vector norm of Matrix.

Wrapper around Matrix::norm - sqrt of sum of squared elements

Template Parameters
MtMatrix type
Parameters
rhsMatrix to compute norm of
Returns
auto

Definition at line 2496 of file matrix.h.

◆ operator*()

template<typename Mt1 , typename Mt2 , std::enable_if_t< Mt1::is_matrix() &&Mt2::is_matrix() &&!std::is_same< Mt1, Mt2 >::value, int > = 0, typename R = hila::type_mul<hila::number_type<Mt1>, hila::number_type<Mt2>>, int n = Mt1::rows(), int m = Mt2::columns()>
Matrix< n, m, R > operator* ( const Mt1 &  a,
const Mt2 &  b 
)
inline

Multiplication operator.

Multiplication type depends on the original types of the multiplied matrices. Defined for the following operations.

Matrix * Matrix:

Standard matrix multiplication where LHS columns must match RHS rows

M.random();
N.random();
auto S = N * M; // Results in a 3 x 3 Matrix since N.rows() = 3 and M.columns = 3
Mtype & random()
dot with matrix - matrix
Definition matrix.h:1296

RowVector * Vector / Vector * RowVector:

Defined as standard dot product between vectors as long as vectors are of same length

//
// Fill V and W
//
auto S = V*W; // Results in MyType scalar

Vector * RowVector is same as outer product which is equivalent to a matrix multiplication

auto S = W*V // Result in n x n Matrix of type MyType

Matrix * Scalar / Scalar * Matrix:

Multiplication operator between Matrix and Scalar are defined in the usual way (element wise)

M.fill(1);
auto S = M*2; // Resulting Matrix is uniformly 2
const auto & fill(const S rhs)
Matrix fill.
Definition matrix.h:937
Template Parameters
Mt1Matrix type for a
Mt2Matrix type for b
nNumber of rows
mNumber of columns
Parameters
aLeft Matrix, Vector or Scalar
bRight Matrix, Vector or Scalar
Returns
Matrix<n, m, R>

Definition at line 2187 of file matrix.h.

◆ operator+()

template<typename Mtype1 , typename Mtype2 , std::enable_if_t< Mtype1::is_matrix() &&Mtype2::is_matrix(), int > = 0, typename Rtype = hila::mat_x_mat_type<Mtype1, Mtype2>, std::enable_if_t<!std::is_same< Mtype1, Rtype >::value, int > = 0>
Rtype operator+ ( const Mtype1 &  a,
const Mtype2 &  b 
)
inline

Addition operator.

Defined for the following operations

Matrix + Matrix:

Addition operator between matrices is defined in the usual way (element wise).

NOTE: Matrices must share same dimensionality.

M.fill(1);
N.fill(1);
S = M + N; // Resulting matrix is uniformly 2

Scalar + Matrix / Matrix + Scalar:

Addition operator between matrix and scalar is defined in the usual way, where the scalar is added to the diagonal elements.

NOTE: Exact definition exist in overloaded functions that can be viewed in source code.

\( M + 2 = M + 2\cdot\mathbb{1} \)

M.fill(0);
S = M + 1; // Resulting matrix is identity matrix
Template Parameters
Mtype1Matrix type for a
Mtype2Matrix type for b
Parameters
aLeft matrix or scalar
bRight matrix or scalar
Returns
Rtype Return matrix of compatible type between Mtype1 and Mtype2

Definition at line 1946 of file matrix.h.

◆ operator-()

template<typename Mtype1 , typename Mtype2 , std::enable_if_t< Mtype1::is_matrix() &&Mtype2::is_matrix(), int > = 0, typename Rtype = hila::mat_x_mat_type<Mtype1, Mtype2>>
Rtype operator- ( const Mtype1 &  a,
const Mtype2 &  b 
)
inline

Subtraction operator.

Defined for the following operations

Matrix - Matrix:

Subtraction operator between matrices is defined in the usual way (element wise).

NOTE: Matrices must share same dimensionality.

M.fill(1);
N.fill(1);
S = M - N; // Resulting matrix is uniformly 0

Scalar - Matrix / Matrix - Scalar:

Subtraction operator between matrix and scalar is defined in the usual way, where the scalar is subtracted from the diagonal elements.

NOTE: Exact definition exist in overloaded functions that can be viewed in source code.

\( M - 2 = M - 2\cdot\mathbb{1} \)

M = 2; // M = 2*I
S = M - 1; // Resulting matrix is identity matrix
Template Parameters
Mtype1Matrix type for a
Mtype2Matrix type for b
Parameters
aLeft matrix or scalar
bRight matrix or scalar
Returns
Rtype Return matrix of compatible type between Mtype1 and Mtype2

Definition at line 2022 of file matrix.h.

◆ operator/()

template<typename Mt , typename S , std::enable_if_t<(Mt::is_matrix() &&hila::is_complex_or_arithmetic< S >::value), int > = 0, typename Rt = hila::mat_scalar_type<Mt, S>>
Rt operator/ ( const Mt &  mat,
const S &  rhs 
)
inline

Division operator.

Defined for following operations

__ Matrix / Scalar: __

Division operator between Matrix and Scalar are defined in the usual way (element wise)

M.fill(2);
auto S = M/2; // Resulting matrix is uniformly 1
Template Parameters
MtMatrix type
SScalar type
Parameters
matMatrix to divide scalar with
rhsScalar to divide with
Returns
Rt Resulting Matrix

Definition at line 2285 of file matrix.h.

◆ operator<<()

template<int n, int m, typename T , typename MT >
std::ostream & operator<< ( std::ostream &  strm,
const Matrix_t< n, m, T, MT > &  A 
)

Stream operator.

Naive Stream operator for formatting Matrix for printing. Simply puts elements one after the other in row major order

Template Parameters
n
m
T
MT
Parameters
strm
A
Returns
std::ostream&

Definition at line 2358 of file matrix.h.

◆ real()

template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto real ( const Mtype arg)
inline

Return real of Matrix or Vector.

Wrapper around Matrix::real

Can be called as:

M_real = real(M);
Array< n, m, hila::arithmetic_type< T > > real(const Array< n, m, T > &arg)
Return real part of Array.
Definition array.h:688
Template Parameters
MtypeMatrix type
Parameters
argObject to compute real part of
Returns
auto Mtype real part of arg

Definition at line 1874 of file matrix.h.

◆ squarenorm()

template<typename Mt , std::enable_if_t< Mt::is_matrix(), int > = 0>
auto squarenorm ( const Mt &  rhs)
inline

Returns square norm of Matrix.

Wrapper around Matrix::squarenorm - sum of squared elements

Can be called as:

auto a = squarenorm(M);
hila::arithmetic_type< T > squarenorm(const Array< n, m, T > &rhs)
Return square norm of Array.
Definition array.h:1018
Template Parameters
MtMatrix type
Parameters
rhsMatrix to compute square norm of
Returns
auto

Definition at line 2482 of file matrix.h.

◆ trace()

template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto trace ( const Mtype arg)
inline

Return trace of square Matrix.

Wrapper around Matrix::trace

Can be called as:

T = trace(M);
Complex definition.
Definition cmplx.h:56
Template Parameters
MtypeMatrix type
Parameters
argObject to compute trace of
Returns
auto Trace of Matrix

Definition at line 1854 of file matrix.h.

◆ transpose()

template<typename Mtype , std::enable_if_t< Mtype::is_matrix(), int > = 0>
auto transpose ( const Mtype arg)
inline

Return transposed Matrix or Vector.

Wrapper around Matrix::transpose

Can be called as:

Matrix<n,m,MyType> M,M_transposed;
M_transposed = transpose(M);
Template Parameters
MtypeMatrix type
Parameters
argObject to be transposed
Returns
auto Mtype transposed matrix

Definition at line 1757 of file matrix.h.