CASM  1.1.0
A Clusters Approach to Statistical Mechanics
LinearAlgebra

Detailed Description

Linear algebra routines.

Functions

void CASM::get_Hermitian (Eigen::MatrixXcd &original_mat, Eigen::MatrixXcd &hermitian_mat, Eigen::MatrixXcd &antihermitian_mat)
 
bool CASM::is_Hermitian (Eigen::MatrixXcd &mat)
 
std::pair< Eigen::MatrixXi, Eigen::MatrixXi > CASM::hermite_normal_form (const Eigen::MatrixXi &M)
 Return the hermite normal form, M == H*V. More...
 
double CASM::angle (const Eigen::Ref< const Eigen::Vector3d > &a, const Eigen::Ref< const Eigen::Vector3d > &b)
 Get angle, in radians, between two vectors on range [0,pi]. More...
 
double CASM::signed_angle (const Eigen::Ref< const Eigen::Vector3d > &a, const Eigen::Ref< const Eigen::Vector3d > &b, const Eigen::Ref< const Eigen::Vector3d > &pos_ref)
 signed angle, in radians, between -pi and pi that describe separation in direction of two vectors More...
 
Eigen::MatrixXd CASM::pretty (const Eigen::MatrixXd &M, double tol)
 Round entries that are within tol of being integer to that integer value. More...
 
template<typename Derived >
Derived::Scalar CASM::triple_product (const Derived &vec0, const Derived &vec1, const Derived &vec2)
 
template<typename Derived >
bool CASM::is_integer (const Eigen::MatrixBase< Derived > &M, double tol)
 Check if Eigen::Matrix is integer. More...
 
template<typename Derived >
bool CASM::is_unimodular (const Eigen::MatrixBase< Derived > &M, double tol)
 Check if Eigen::Matrix is unimodular. More...
 
template<typename Derived >
bool CASM::is_diagonal (const Eigen::MatrixBase< Derived > &M, double tol=TOL)
 Check if Eigen::Matrix is diagonal. More...
 
template<typename Derived >
Eigen::CwiseUnaryOp< decltype(Local::round_l< typename Derived::Scalar >), const Derived > CASM::round (const Eigen::MatrixBase< Derived > &val)
 Round Eigen::MatrixXd. More...
 
template<typename Derived >
Eigen::CwiseUnaryOp< decltype(Local::iround_l< typename Derived::Scalar >), const Derived > CASM::iround (const Eigen::MatrixBase< Derived > &val)
 Round Eigen::MatrixXd to Eigen::MatrixXi. More...
 
template<typename Derived >
Eigen::CwiseUnaryOp< decltype(Local::lround_l< typename Derived::Scalar >), const Derived > CASM::lround (const Eigen::MatrixBase< Derived > &val)
 Round Eigen::MatrixXd to Eigen::MatrixXl. More...
 
template<typename Derived >
Derived::Scalar CASM::matrix_minor (const Eigen::MatrixBase< Derived > &M, int row, int col)
 Return the minor of integer Matrix M element row, col. More...
 
template<typename Derived >
Eigen::Matrix< typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime > CASM::cofactor (const Eigen::MatrixBase< Derived > &M)
 Return cofactor matrix. More...
 
template<typename Derived >
Eigen::Matrix< typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime > CASM::adjugate (const Eigen::MatrixBase< Derived > &M)
 Return adjugate matrix. More...
 
template<typename Derived >
Eigen::Matrix< typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime > CASM::inverse (const Eigen::MatrixBase< Derived > &M)
 Return the integer inverse matrix of an invertible integer matrix. More...
 
template<typename DerivedIn , typename DerivedOut >
void CASM::smith_normal_form (const Eigen::MatrixBase< DerivedIn > &M, Eigen::MatrixBase< DerivedOut > &U, Eigen::MatrixBase< DerivedOut > &S, Eigen::MatrixBase< DerivedOut > &V)
 Return the smith normal form, M == U*S*V. More...
 
Eigen::Matrix3d CASM::polar_decomposition (Eigen::Matrix3d const &mat)
 
std::vector< Eigen::Matrix3i > CASM::_unimodular_matrices (bool positive, bool negative, int range=1)
 
const std::vector< Eigen::Matrix3i > & CASM::positive_unimodular_matrices ()
 
const std::vector< Eigen::Matrix3i > & CASM::negative_unimodular_matrices ()
 
template<int range = 1>
const std::vector< Eigen::Matrix3i > & CASM::unimodular_matrices ()
 

Function Documentation

◆ _unimodular_matrices()

std::vector< Eigen::Matrix3i > CASM::_unimodular_matrices ( bool  positive,
bool  negative,
int  range = 1 
)

Definition at line 646 of file CASM_Eigen_math.cc.

◆ adjugate()

template<typename Derived >
Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime> CASM::adjugate ( const Eigen::MatrixBase< Derived > &  M)

Return adjugate matrix.

Returns
adjugate matrix
Parameters
Mmatrix

The adjugate matrix is the transpose of the cofactor matrix. The adjugate matrix is related to the inverse of a matrix through

M.inverse() == adjugate(M)/M.determinant()
Eigen::Matrix< typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime > adjugate(const Eigen::MatrixBase< Derived > &M)
Return adjugate matrix.

Definition at line 320 of file CASM_Eigen_math.hh.

◆ angle()

double CASM::angle ( const Eigen::Ref< const Eigen::Vector3d > &  a,
const Eigen::Ref< const Eigen::Vector3d > &  b 
)

Get angle, in radians, between two vectors on range [0,pi].

Definition at line 611 of file CASM_Eigen_math.cc.

◆ cofactor()

template<typename Derived >
Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime> CASM::cofactor ( const Eigen::MatrixBase< Derived > &  M)

Return cofactor matrix.

Returns
cofactor matrix
Parameters
Mmatrix

The cofactor matrix is

C(i,j) = pow(-1,i+j)*matrix_minor(M,i,j)
Derived::Scalar matrix_minor(const Eigen::MatrixBase< Derived > &M, int row, int col)
Return the minor of integer Matrix M element row, col.

Definition at line 292 of file CASM_Eigen_math.hh.

◆ get_Hermitian()

void CASM::get_Hermitian ( Eigen::MatrixXcd &  original_mat,
Eigen::MatrixXcd &  hermitian_mat,
Eigen::MatrixXcd &  antihermitian_mat 
)

Given a matrix, original_mat, this function calculates both the Hermitian and anti-Hermitian parts.

Parameters
[in]original_matThe matrix we want to get the Hermitian and anti-Hermitian parts of
[in,out]hermitian_matThe Hermitian matrix
[in,out]antihermitian_matThe anti-Hermitian matrix

Definition at line 486 of file CASM_Eigen_math.cc.

◆ hermite_normal_form()

std::pair< Eigen::MatrixXi, Eigen::MatrixXi > CASM::hermite_normal_form ( const Eigen::MatrixXi &  M)

Return the hermite normal form, M == H*V.

Parameters
Ma square integer matrix
Exceptions
std::runtime_errorif M is not full rank
Returns
std::pair<H,V>

H set to upper triangular matrix H with:

  • H(i,i) > 0
  • 0 <= H(r,i) < H(i,i) for r > i
  • H(r,i) == 0 for r < i V set to unimodular matrix V
  • V.determinant() = +/- 1

Definition at line 524 of file CASM_Eigen_math.cc.

◆ inverse()

template<typename Derived >
Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime> CASM::inverse ( const Eigen::MatrixBase< Derived > &  M)

Return the integer inverse matrix of an invertible integer matrix.

Returns
a 3x3 integer matrix inverse
Parameters
Man 3x3 invertible integer matrix

Definition at line 359 of file CASM_Eigen_math.hh.

◆ iround()

template<typename Derived >
Eigen::CwiseUnaryOp<decltype(Local::iround_l<typename Derived::Scalar>), const Derived> CASM::iround ( const Eigen::MatrixBase< Derived > &  val)

Round Eigen::MatrixXd to Eigen::MatrixXi.

Returns
an Eigen:MatrixXi
Parameters
MEigen::MatrixXd to be rounded to integer

For each coefficient, sets

Mint(i,j) = boost::math::iround(Mdouble(i,
j))
Eigen::CwiseUnaryOp< decltype(Local::iround_l< typename Derived::Scalar >), const Derived > iround(const Eigen::MatrixBase< Derived > &val)
Round Eigen::MatrixXd to Eigen::MatrixXi.

Definition at line 218 of file CASM_Eigen_math.hh.

◆ is_diagonal()

template<typename Derived >
bool CASM::is_diagonal ( const Eigen::MatrixBase< Derived > &  M,
double  tol = TOL 
)

Check if Eigen::Matrix is diagonal.

Definition at line 186 of file CASM_Eigen_math.hh.

◆ is_Hermitian()

bool CASM::is_Hermitian ( Eigen::MatrixXcd &  mat)

Returns true if the specified matrix is Hermitian.

Checks if the matrix is the same as the conjugate

Definition at line 503 of file CASM_Eigen_math.cc.

◆ is_integer()

template<typename Derived >
bool CASM::is_integer ( const Eigen::MatrixBase< Derived > &  M,
double  tol 
)

Check if Eigen::Matrix is integer.

Definition at line 166 of file CASM_Eigen_math.hh.

◆ is_unimodular()

template<typename Derived >
bool CASM::is_unimodular ( const Eigen::MatrixBase< Derived > &  M,
double  tol 
)

Check if Eigen::Matrix is unimodular.

Definition at line 178 of file CASM_Eigen_math.hh.

◆ lround()

template<typename Derived >
Eigen::CwiseUnaryOp<decltype(Local::lround_l<typename Derived::Scalar>), const Derived> CASM::lround ( const Eigen::MatrixBase< Derived > &  val)

Round Eigen::MatrixXd to Eigen::MatrixXl.

Returns
an Eigen:MatrixXl
Parameters
MEigen::MatrixXd to be rounded to integer

For each coefficient, sets

Mint(i,j) = std::lround(Mdouble(i, j))
Eigen::CwiseUnaryOp< decltype(Local::lround_l< typename Derived::Scalar >), const Derived > lround(const Eigen::MatrixBase< Derived > &val)
Round Eigen::MatrixXd to Eigen::MatrixXl.

Definition at line 234 of file CASM_Eigen_math.hh.

◆ matrix_minor()

template<typename Derived >
Derived::Scalar CASM::matrix_minor ( const Eigen::MatrixBase< Derived > &  M,
int  row,
int  col 
)

Return the minor of integer Matrix M element row, col.

Returns
the minor of element row, col
Parameters
Mmatrix
row,colelement row and column

The minor of element row, col of a matrix is the determinant of the submatrix of M which does not include any elements from row 'row' or column 'col'

Definition at line 250 of file CASM_Eigen_math.hh.

◆ negative_unimodular_matrices()

const std::vector< Eigen::Matrix3i > & CASM::negative_unimodular_matrices ( )

Definition at line 680 of file CASM_Eigen_math.cc.

◆ polar_decomposition()

Eigen::Matrix3d CASM::polar_decomposition ( Eigen::Matrix3d const &  mat)

Definition at line 641 of file CASM_Eigen_math.cc.

◆ positive_unimodular_matrices()

const std::vector< Eigen::Matrix3i > & CASM::positive_unimodular_matrices ( )

Definition at line 674 of file CASM_Eigen_math.cc.

◆ pretty()

Eigen::MatrixXd CASM::pretty ( const Eigen::MatrixXd &  M,
double  tol 
)

Round entries that are within tol of being integer to that integer value.

Definition at line 629 of file CASM_Eigen_math.cc.

◆ round()

template<typename Derived >
Eigen::CwiseUnaryOp<decltype(Local::round_l<typename Derived::Scalar>), const Derived> CASM::round ( const Eigen::MatrixBase< Derived > &  val)

Round Eigen::MatrixXd.

Returns
an Eigen:MatrixXd
Parameters
MEigen::MatrixXd to be rounded

For each coefficient, sets

M(i,j) = boost::math::round(Mdouble(i, j))
Eigen::CwiseUnaryOp< decltype(Local::round_l< typename Derived::Scalar >), const Derived > round(const Eigen::MatrixBase< Derived > &val)
Round Eigen::MatrixXd.

Definition at line 202 of file CASM_Eigen_math.hh.

◆ signed_angle()

double CASM::signed_angle ( const Eigen::Ref< const Eigen::Vector3d > &  a,
const Eigen::Ref< const Eigen::Vector3d > &  b,
const Eigen::Ref< const Eigen::Vector3d > &  pos_ref 
)

signed angle, in radians, between -pi and pi that describe separation in direction of two vectors

return signed angle, in radians, between -pi and pi that describe separation in direction of two vectors

Definition at line 618 of file CASM_Eigen_math.cc.

◆ smith_normal_form()

template<typename DerivedIn , typename DerivedOut >
void CASM::smith_normal_form ( const Eigen::MatrixBase< DerivedIn > &  M,
Eigen::MatrixBase< DerivedOut > &  U,
Eigen::MatrixBase< DerivedOut > &  S,
Eigen::MatrixBase< DerivedOut > &  V 
)

Return the smith normal form, M == U*S*V.

Parameters
M3x3 integer matrix
[out]Uset to U: is integer, U.determinant() == 1
[out]Sset to S: is diagonal, integer, and nonnegative, S(i,i) divides S(i+1,i+1) for all i
[out]Vset to V: is integer, V.determinant() == 1

Adapted from Matlab implementation written by John Gilbert (gilbe.nosp@m.rt@p.nosp@m.arc.x.nosp@m.erox.nosp@m..com):

Definition at line 376 of file CASM_Eigen_math.hh.

◆ triple_product()

template<typename Derived >
Derived::Scalar CASM::triple_product ( const Derived &  vec0,
const Derived &  vec1,
const Derived &  vec2 
)

Definition at line 158 of file CASM_Eigen_math.hh.

◆ unimodular_matrices()

template<int range = 1>
const std::vector<Eigen::Matrix3i>& CASM::unimodular_matrices ( )

Definition at line 528 of file CASM_Eigen_math.hh.