CASM  1.1.0
A Clusters Approach to Statistical Mechanics
DoFSet.cc
Go to the documentation of this file.
3 #include <unordered_set>
4 #include <vector>
5 
7 
8 namespace CASM {
9 namespace xtal {
12 std::vector<std::string> component_descriptions(DoFSet const &dofset) {
13  if (dofset.basis().isIdentity(TOL) &&
14  dofset.basis().rows() == dofset.basis().cols())
15  return dofset.traits().variable_descriptions();
16  else
17  return dofset.component_names();
18 }
19 
20 bool DoFSetIsEquivalent_f::_traits_match(const DoFSet &other_value) const {
21  return m_reference_dofset.traits().name() == other_value.traits().name();
22 }
23 
25  const DoFSet &other_value) const {
26  const auto &reference_basis = m_reference_dofset.basis();
27  const auto &other_basis = other_value.basis();
28 
29  if (reference_basis.cols() != other_basis.cols()) return false;
30 
31  // If not a square matrix, make sure that column-space _after_basis is
32  // identical to that of _before_basis
33  if (reference_basis.cols() < other_basis.rows()) {
34  // Find rank of augmented matrix. If it is the same as rank of DoF Basis,
35  // then the two matrices are similar (and, thus, equivalent)
36  Eigen::MatrixXd aug(reference_basis.rows(), 2 * reference_basis.cols());
37  aug << reference_basis, other_basis;
38 
39  Eigen::ColPivHouseholderQR<Eigen::MatrixXd> qr(aug);
40  qr.setThreshold(m_tol);
41  if (qr.rank() != m_reference_dofset.dim()) return false;
42  }
43 
44  return true;
45 }
46 
47 bool DoFSetIsEquivalent_f::operator()(const DoFSet &other_value) const {
48  return this->_traits_match(other_value) &&
49  this->_basis_spans_same_space(other_value);
50 }
51 
53  const Eigen::MatrixXd &to_basis,
54  double tol) {
55  Eigen::MatrixXd U = from_basis.colPivHouseholderQr().solve(to_basis);
56  if (!(U.transpose() * U).eval().isIdentity(tol)) {
57  throw std::runtime_error("Cannot find orthogonal symmetry representation!");
58  }
59  return U;
60 }
61 } // namespace xtal
62 } // namespace CASM
63 
64 namespace CASM {
65 namespace sym {
68  Eigen::MatrixXd transformation_matrix = dof.traits().symop_to_matrix(
70  return xtal::DoFSet(dof.traits(), dof.component_names(),
71  transformation_matrix * dof.basis());
72 }
73 } // namespace sym
74 
75 } // namespace CASM
std::string const & name() const
Const access of name.
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.
std::vector< std::string > const & variable_descriptions() const
Returns expanded description of each standard_var_name.
Eigen::MatrixXd const & basis() const
Matrix that relates DoFSet variables to a conventional coordiante system.
Definition: DoFSet.hh:67
Index dim() const
Definition: DoFSet.hh:64
BasicTraits const & traits() const
Returns traits object for the DoF type of this DoFSet.
Definition: DoFSet.hh:60
const std::vector< std::string > & component_names() const
Returns the names of each of the component axes.
Definition: DoFSet.hh:55
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
xtal::Coordinate copy_apply(const xtal::SymOp &op, xtal::Coordinate coord)
Copy and apply SymOp to a Coordinate.
Definition: Coordinate.cc:354
std::vector< std::string > component_descriptions(DoFSet const &dofset)
Definition: DoFSet.cc:12
Eigen::MatrixXd dofset_transformation_matrix(const Eigen::MatrixXd &from_basis, const Eigen::MatrixXd &to_basis, double tol)
Create the symmtery representation for going from one basis to another.
Definition: DoFSet.cc:52
Main CASM namespace.
Definition: APICommand.hh:8
Eigen::MatrixXd MatrixXd
const double TOL
Definition: definitions.hh:30
double m_tol
Tolerance value for making comparisons.
Definition: DoFSet.hh:163
bool _traits_match(const DoFSet &other_value) const
Returns true if the traits match. Only compares the names.
Definition: DoFSet.cc:20
bool operator()(const DoFSet &other_value) const
Definition: DoFSet.cc:47
bool _basis_spans_same_space(const DoFSet &other_value) const
Definition: DoFSet.cc:24
DoFSet m_reference_dofset
Values passed to operator() will be compared against this.
Definition: DoFSet.hh:160