CASM  1.1.0
A Clusters Approach to Statistical Mechanics
DoFIsEquivalent.cc
Go to the documentation of this file.
2 
3 #include "casm/basis_set/DoF.hh"
6 
7 namespace CASM {
8 
9 bool DoFIsEquivalent::operator()(DoFSet const &other) const {
10  return _label_equiv(other) && _vector_equiv(other.basis(), m_dof.basis());
11 }
12 
13 bool DoFIsEquivalent::operator()(xtal::SymOp const &_op) const {
14  return _vector_equiv(
15  m_dof.basis(),
18  get_time_reversal(_op)) *
19  m_dof.basis());
20 }
21 
23  DoFSet const &other) const {
24  return _label_equiv(other) &&
26  other.basis(),
29  get_time_reversal(_op)) *
30  m_dof.basis());
31 }
32 
33 bool DoFIsEquivalent::_label_equiv(DoFSet const &other) const {
34  return m_dof.type_name() == m_dof.type_name() && m_dof.size() == other.size();
35 }
36 
38  Eigen::Ref<const Eigen::MatrixXd> const &_before_basis,
39  Eigen::Ref<const Eigen::MatrixXd> const &_after_basis) const {
40  if (_before_basis.cols() != _after_basis.cols()) return false;
41 
42  // If not a square matrix, make sure that column-space _after_basis is
43  // identical to that of _before_basis
44  if (_before_basis.cols() < _before_basis.rows()) {
45  // Find rank of augmented matrix. If it is the same as rank of DoF Basis,
46  // then the two matrices are similar (and, thus, equivalent)
47  Eigen::MatrixXd aug(_before_basis.rows(), 2 * _before_basis.cols());
48  aug << _before_basis, _after_basis;
49  // std::cout << "Basis before: \n" << _before_basis << "\n\n"
50  // << "Basis after: \n" << _after_basis << "\n\n"
51  // << "DoF size: " << m_dof.size() << "\n\n"
52  // << "Augmented rank: " << aug.colPivHouseholderQr().rank() <<
53  // "\n\n";
54  Eigen::ColPivHouseholderQR<Eigen::MatrixXd> qr(aug);
55  qr.setThreshold(m_tol);
56  if (qr.rank() != m_dof.size()) return false;
57  }
58 
59  // Find the transform m_U where
60  // _before_basis*m_U = _after_basis
61  m_U = _before_basis.colPivHouseholderQr().solve(_after_basis);
62  if (!(m_U.transpose() * m_U).eval().isIdentity(m_tol)) {
63  throw std::runtime_error(
64  "Cannot find orthogonal symmetry representation for DoF \"" +
65  m_dof.type_name() + "\". Please review inputs.\n");
66  }
67  return true;
68 }
69 
70 } // namespace CASM
Specifies traits of (possibly) anisotropic crystal properties.
Eigen::MatrixXd symop_to_matrix(Eigen::Ref< const Eigen::Matrix3d > const &_matrix, Eigen::Ref< const Eigen::Vector3d > const &_tau, bool time_reversal) const
Generate a symmetry representation for the supporting vector space.
bool _label_equiv(DoFSet const &_other) const
returns true if m_dof and _other are same type and same size
bool operator()(DoFSet const &_other) const
bool _vector_equiv(Eigen::Ref< const Eigen::MatrixXd > const &_before_basis, Eigen::Ref< const Eigen::MatrixXd > const &_after_basis) const
Index size() const
Returns the number of components in this DoFSet.
Definition: DoFSet.hh:191
Eigen::MatrixXd const & basis() const
Definition: DoFSet.hh:268
std::string const & type_name() const
Definition: DoFSet.hh:197
const SymOp::matrix_type & get_matrix(const SymOp &op)
Definition: SymOp.cc:167
bool get_time_reversal(const SymOp &op)
Definition: SymOp.cc:171
const SymOp::vector_type & get_translation(const SymOp &op)
Definition: SymOp.cc:169
Main CASM namespace.
Definition: APICommand.hh:8
Eigen::MatrixXd MatrixXd