HILA
|
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"
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 Mt1 , typename Mt2 , typename Mt3 , std::enable_if_t< Mt1::is_matrix() &&Mt2::is_matrix() &&Mt3::is_matrix(), int > = 0> | |
void | mult (const Mt1 &a, const Mt2 &b, Mt3 &res) |
compute product of two matrices and write result to existing matrix | |
template<typename Mt1 , typename Mt2 , typename Mt3 , std::enable_if_t< Mt1::is_matrix() &&Mt2::is_matrix() &&Mt3::is_matrix(), int > = 0> | |
void | mult_aa (const Mt1 &a, const Mt2 &b, Mt3 &res) |
compute hermitian conjugate of product of two matrices and write result to existing matrix | |
template<typename Mt1 , typename S , typename Mt2 , std::enable_if_t<(Mt1::is_matrix() &&Mt2::is_matrix() &&hila::is_complex_or_arithmetic< S >::value), int > = 0> | |
void | mult (const Mt1 &a, const S &b, Mt2 &res) |
compute product of a matrix and a scalar and write result to existing matrix | |
template<typename S , typename Mt1 , typename Mt2 , std::enable_if_t<(Mt1::is_matrix() &&Mt2::is_matrix() &&hila::is_complex_or_arithmetic< S >::value), int > = 0> | |
void | mult (const S &a, const Mt1 &b, Mt2 &res) |
compute product of a scalar and a matrix and write result to existing matrix | |
template<typename Mt1 , typename Mt2 , typename Mt3 , std::enable_if_t< Mt1::is_matrix() &&Mt2::is_matrix() &&Mt3::is_matrix(), int > = 0> | |
void | mult_add (const Mt1 &a, const Mt2 &b, Mt3 &res) |
compute product of two matrices and add result to existing matrix | |
template<typename Mt1 , typename S , typename Mt2 , std::enable_if_t<(Mt1::is_matrix() &&Mt2::is_matrix() &&hila::is_complex_or_arithmetic< S >::value), int > = 0> | |
void | mult_add (const Mt1 &a, const S &b, Mt2 &res) |
compute product of a matrix and a scalar and add result to existing matrix | |
template<typename S , typename Mt1 , typename Mt2 , std::enable_if_t<(Mt1::is_matrix() &&Mt2::is_matrix() &&hila::is_complex_or_arithmetic< S >::value), int > = 0> | |
void | mult_add (const S &a, const Mt1 &b, Mt2 &res) |
compute product of a scalar and a matrix and add result to existing matrix | |
template<typename Mt1 , typename Mt2 , typename Mt3 , std::enable_if_t< Mt1::is_matrix() &&Mt2::is_matrix() &&Mt3::is_matrix(), int > = 0> | |
void | mult_sub (const Mt1 &a, const Mt2 &b, Mt3 &res) |
compute product of two matrices and subtract result from existing matrix | |
template<typename Mt1 , typename S , typename Mt2 , std::enable_if_t<(Mt1::is_matrix() &&Mt2::is_matrix() &&hila::is_complex_or_arithmetic< S >::value), int > = 0> | |
void | mult_sub (const Mt1 &a, const S &b, Mt2 &res) |
compute product of a matrix and a scalar and subtract result from existing matrix | |
template<typename S , typename Mt1 , typename Mt2 , std::enable_if_t<(Mt1::is_matrix() &&Mt2::is_matrix() &&hila::is_complex_or_arithmetic< S >::value), int > = 0> | |
void | mult_sub (const S &a, const Mt1 &b, Mt2 &res) |
compute product of a scalar and a matrix and subtract result from 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 , typename atype = hila::arithmetic_type<T>> | |
Matrix_t< n, m, T, MT > | altexp (const Matrix_t< n, m, T, MT > &mat, int &niter) |
Calculate exp of a square matrix. | |
template<int n, int m, typename T , typename MT , typename atype = hila::arithmetic_type<T>> | |
Matrix_t< n, m, T, MT > | altexp (const Matrix_t< n, m, T, MT > &mat) |
Calculate exp of a square matrix. | |
template<int n, int m, typename T , typename MT > | |
void | mult_exp (const Matrix_t< n, m, T, MT > &mat, const Matrix_t< n, m, T, MT > &mmat, Matrix_t< n, m, T, MT > &r, Matrix_t< n, m, T, MT > &dr, const int order=20) |
Calculate mmat*exp(mat) and trace(mmat*dexp(mat)) | |
template<int n, int m, typename T , typename MT , typename Mt , std::enable_if_t< Mt::is_matrix() &&Mt::is_square() &&Mt::rows()==n, int > = 0> | |
int | chexp (const Matrix_t< n, m, T, MT > &mat, Matrix_t< n, m, T, MT > &omat, Mt(&pl)[n]) |
Calculate exp of a square matrix. | |
template<int n, int m, typename T , typename MT , typename Mt , std::enable_if_t< Mt::is_matrix() &&Mt::is_square() &&Mt::rows()==n, int > = 0> | |
void | chexp (const Matrix_t< n, m, T, MT > &mat, Matrix_t< n, m, T, MT > &omat, Mt(&domat)[n][m]) |
Calculate exp and dexp of a square matrix. | |
template<int n, int m, typename T , typename MT > | |
void | mult_chexp (const Matrix_t< n, m, T, MT > &mat, const Matrix_t< n, m, T, MT > &mmat, Matrix_t< n, m, T, MT > &omat, Matrix_t< n, m, T, MT > &domat) |
Calculate exp(mat).dagger()*mmat*exp(mat) and trace(exp(mat).dagger*mmat*dexp(mat)) | |
template<int n, int m, typename T , typename MT > | |
void | chexpk (const Matrix_t< n, m, T, MT > &mat, Matrix_t< n, m, T, MT > &omat, Matrix_t< n, m, T, MT > &kmats) |
Calculate exp(mat) and the decomposition k_{i,j} of dexp in terms bilinears of powers of mat. | |
template<int n, int m, typename T , typename MT > | |
void | mult_chexpk_fast (const Matrix_t< n, m, T, MT > &mat, const Matrix_t< n, m, T, MT > &texp, const Matrix_t< n, m, T, MT > &kmats, const Matrix_t< n, m, T, MT > &mmat, Matrix_t< n, m, T, MT > &omat, Matrix_t< n, m, T, MT > &domat) |
Calculate exp(mat).dagger()*mmat*exp(mat) and trace(exp(mat).dagger*mmat*dexp(mat)) from output of chexpk. | |
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. | |
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.
|
inline |
Return absolute value Matrix or Vector.
Wrapper around Matrix::abs
Can be called as:
Mtype | Matrix type |
arg | Object to compute absolute value of |
|
inline |
Return adjoint Matrix.
Wrapper around Matrix::adjoint
Can be called as:
Mtype | Matrix type |
arg | Object to compute adjoint of |
|
inline |
Calculate exp of a square matrix.
Adds hihger order terms till matrix norm converges within floating point accuracy
mat | Matrix to compute exponential for |
|
inline |
Calculate exp of a square matrix.
Adds higher order terms till matrix norm converges within floating point accuracy and returns required number of iterations
mat | Matrix to compute exponential for |
niter | (output) number of iteration performed till converges was reached |
|
inline |
Calculate exp and dexp of a square matrix.
Computation is done using iterative Cayley-Hamilton (cf. from arXiv:2404.07704)
mat | Matrix to compute exponential for |
omat | Matrix to which exponential of mat gets stored |
domat | matrix of Matrices to which the derivatives of the exponential with respect to the components of mat gets stored |
|
inline |
Calculate exp of a square matrix.
Computation is done using iterative Cayley-Hamilton (cf. from arXiv:2404.07704)
mat | Matrix to compute exponential for |
omat | Matrix to which exponential of mat gets stored (optional) |
pl | array of n+1 temporary nxn Matrices (optional) |
|
inline |
Calculate exp(mat) and the decomposition k_{i,j} of dexp in terms bilinears of powers of mat.
Computation is done using iterative Cayley-Hamilton (cf. from arXiv:2404.07704) exp(mat)^a_b = \sum_{i=0}^{n-1} r_{i} (mat^i)^a_b dexp(X)^a_b/dX^c_d|_X=mat = \sum_{i,j=0}^{n-1} k_{i,j} (mat^i)^d_b (mat^j)^a_c
n | Number of rows of Matrix |
m | Number of columns of Matrix |
T | Matrix element type |
MT | Matrix type |
mat | Matrix of which to compute exponential |
omat | Matrix to which exponential of mat gets stored |
kmats | Matrix to which decomposition coefficients of dexp/dX|X=mat in terms of bilinears of powers of mat get stroed |
|
inline |
Calculate exp of a square matrix.
Computation is done using iterative Cayley-Hamilton (cf. from arXiv:2404.07704) with minimal temporary storage
mat | Matrix to compute exponential for |
|
inline |
Return conjugate Matrix or Vector.
Wrapper around Matrix::conj
Can be called as:
Mtype | Matrix type |
arg | Object to be conjugated |
|
inline |
Return dagger of Matrix.
Wrapper around Matrix::adjoint
Same as adjoint
Mtype | Matrix type |
arg | Object to compute dagger of |
|
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
mat | Matrix to compute exponential for |
order | order to compute exponential to |
|
inline |
Return imaginary of Matrix or Vector.
Wrapper around Matrix::imag
Can be called as:
Mtype | Matrix type |
arg | Object to compute imag part of |
|
inline |
Returns trace of product of two matrices.
Wrapper around Matrix::mul_trace in the form
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Calculate exp(mat).dagger()*mmat*exp(mat) and trace(exp(mat).dagger*mmat*dexp(mat))
exp and dexp computation is done using iterative Cayley-Hamilton (arXiv:2404.07704) Calculate omat[i][j] = (exp(mat).dagger() * mmat * exp(mat))[i][j] and domat[i][j] = trace(exp(mat).dagger() * mmat * dexp(mat)/dmat[j][i]) for given matrices mat and mmat, using iterative Cayley-Hamilton method (arXiv:2404.07704) for computing exp(mat) and dexp(mat)/dmat[][]
mat | Matrix to compute exponential for |
mmat | Matrix to multiply with |
omat | Matrix to which exp(mat).dagger()*mmat*exp(mat) gets stored |
domat | matrix to which trace(exp(mat).dagger()*mmat*dexp(mat)) gets stored |
|
inline |
Calculate exp(mat).dagger()*mmat*exp(mat) and trace(exp(mat).dagger*mmat*dexp(mat)) from output of chexpk.
exp is given and dexp computation is done using decompisition matrix kmats_{i,j}, as computed by chexpk using iterative Cayley-Hamilton (arXiv:2404.07704) Calculate omat[i][j] = (exp(mat).dagger() * mmat * exp(mat))[i][j] and domat[i][j] = trace(exp(mat).dagger() * mmat * dexp(mat)/dmat[j][i]) for given matrices mat, exp(mat), kmats, and mmat
mat | Matrix to compute exponential for |
texp | Matrix containing exp(mat) |
kmats | Matrix containing decomposition coefficients k_{i,j} of dexp(X)^a_b/dX^c_d = k_{i,j} (X^i)^d_b (X^j)^a_c |
mmat | Matrix to multiply with |
omat | Matrix to which exp(mat).dagger()*mmat*exp(mat) gets stored |
domat | matrix to which trace(exp(mat).dagger()*mmat*dexp(mat)) gets stored |
|
inline |
Calculate mmat*exp(mat) and trace(mmat*dexp(mat))
exp and dexp computation is done using factorized, truncated Taylor series method
mat | Matrix to compute exponential for |
mmat | Matrix to multiply with |
r | Matrix to which mmat*exp(mat) gets stored |
dr | matrix to which trace(mmat*dexp(mat)) gets stored |
order | = 20 integer input parameter giving the Taylor series truncation order |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Multiplication operator.
Multiplication type depends on the original types of the multiplied matrices. Defined for the following operations.
Standard matrix multiplication where LHS columns must match RHS rows
RowVector * Vector / Vector * RowVector:
Defined as standard dot product between vectors as long as vectors are of same length
Vector * RowVector is same as outer product which is equivalent to a matrix multiplication
Matrix * Scalar / Scalar * Matrix:
Multiplication operator between Matrix and Scalar are defined in the usual way (element wise)
|
inline |
Addition operator.
Defined for the following operations
Addition operator between matrices is defined in the usual way (element wise).
NOTE: Matrices must share same dimensionality.
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} \)
a | Left matrix or scalar |
b | Right matrix or scalar |
|
inline |
Subtraction operator.
Defined for the following operations
Subtraction operator between matrices is defined in the usual way (element wise).
NOTE: Matrices must share same dimensionality.
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} \)
a | Left matrix or scalar |
b | Right matrix or scalar |
|
inline |
Division operator.
Defined for following operations
__ Matrix / Scalar: __
Division operator between Matrix and Scalar are defined in the usual way (element wise)
Mt | Matrix type |
S | Scalar type |
mat | Matrix to divide scalar with |
rhs | Scalar to divide with |
std::ostream & operator<< | ( | std::ostream & | strm, |
const Matrix_t< n, m, T, MT > & | A | ||
) |
|
inline |
Return real of Matrix or Vector.
Wrapper around Matrix::real
Can be called as:
Mtype | Matrix type |
arg | Object to compute real part of |
|
inline |
Returns square norm of Matrix.
Wrapper around Matrix::squarenorm - sum of squared elements
Can be called as:
Mt | Matrix type |
rhs | Matrix to compute square norm of |
|
inline |
Return trace of square Matrix.
Wrapper around Matrix::trace
Can be called as:
Mtype | Matrix type |
arg | Object to compute trace of |
|
inline |
Return transposed Matrix or Vector.
Wrapper around Matrix::transpose
Can be called as:
Mtype | Matrix type |
arg | Object to be transposed |