CASM  1.1.0
A Clusters Approach to Statistical Mechanics
Niggli.hh
Go to the documentation of this file.
1 #ifndef NIGGLI_HH
2 #define NIGGLI_HH
3 
4 #include <vector>
5 
6 #include "casm/external/Eigen/Core"
8 
9 namespace CASM {
10 
11 namespace xtal {
12 class Lattice;
13 
52 class NiggliRep {
53  public:
56  static std::vector<Eigen::Matrix3d> const &cell_invariant_transforms();
57 
58  NiggliRep(const Lattice &init_lat);
59 
60  NiggliRep(const Eigen::Matrix3d &init_lat_col_mat);
61 
63  double A() const;
64 
66  double B() const;
67 
69  double C() const;
70 
72  double ksi() const;
73 
75  double eta() const;
76 
78  double zeta() const;
79 
80  const Eigen::Matrix3d &metrical_matrix() const;
81 
83  bool meets_criteria_1(double compare_tol) const;
84 
86  bool meets_criteria_2(double compare_tol) const;
87 
89  bool meets_criteria_3(double compare_tol) const;
90 
92  bool meets_criteria_4(double compare_tol) const;
93 
95  bool meets_criteria_5(double compare_tol) const;
96 
98  bool meets_criteria_6(double compare_tol) const;
99 
101  bool meets_criteria_7(double compare_tol) const;
102 
105  bool meets_criteria_8(double compare_tol) const;
106 
108  bool is_niggli(double compare_tol) const;
109 
111  bool is_niggli_type1(double compare_tol) const;
112 
114  bool is_niggli_type2(double compare_tol) const;
115 
116  void debug_criteria(double compare_tol) const;
117 
119  Index niggli_index(double compare_tol) const;
120 
121  private:
124 
126  const double m_scale_factor;
127 };
128 
130 Index niggli_index(const Eigen::Matrix3d &test_lat_mat, double compare_tol);
131 
134 Lattice niggli(const Lattice &in_lat, double compare_tol,
135  bool keep_handedness = false);
136 
139 bool is_niggli(const Eigen::Matrix3d &test_lat_mat, double compare_tol);
140 
142 bool is_niggli(const Lattice &test_lat, double compare_tol);
143 
147  double compare_tol);
148 
153  const Eigen::Matrix3d &low_score_lat_mat,
154  Eigen::Matrix3d &high_score_lat_mat, double compare_tol);
155 
158 bool standard_orientation_compare(const Eigen::Matrix3d &low_score_lat_mat,
159  const Eigen::Matrix3d &high_score_lat_mat,
160  double compare_tol);
161 
163  public:
164  StandardOrientationCompare(double _tol) : m_tol(_tol) {}
165 
166  bool operator()(const Eigen::Matrix3d &LHS,
167  const Eigen::Matrix3d &RHS) const {
168  return standard_orientation_compare(LHS, RHS, m_tol);
169  }
170 
171  private:
172  double m_tol;
173 };
175 } // namespace xtal
176 } // namespace CASM
177 
178 #endif
double zeta() const
2ab*cos(gamma)
Definition: Niggli.cc:67
double ksi() const
2bc*cos(alpha)
Definition: Niggli.cc:63
const Eigen::Matrix3d & metrical_matrix() const
Definition: Niggli.cc:69
NiggliRep(const Lattice &init_lat)
Definition: Niggli.cc:54
bool meets_criteria_3(double compare_tol) const
For type I: ksi>0 && eta>0 && zeta>0 (all angles < 90)
Definition: Niggli.cc:108
bool is_niggli(double compare_tol) const
True if all conditions are true, and either 4 OR 3 is false.
Definition: Niggli.cc:232
bool meets_criteria_7(double compare_tol) const
|zeta|<=A OR (zeta==A, eta<=2*ksi) OR (zeta==-A, eta==0)
Definition: Niggli.cc:180
bool meets_criteria_5(double compare_tol) const
|ksi|<=B OR (ksi==B, zeta<=2*eta) OR (ksi==-B, zeta==0)
Definition: Niggli.cc:133
double eta() const
2ac*cos(beta)
Definition: Niggli.cc:65
Index niggli_index(double compare_tol) const
Number of criteria met.
Definition: Niggli.cc:257
bool meets_criteria_4(double compare_tol) const
For type II: ksi<=0 && eta<=0 && zeta<=0 (all angles >= 90)
Definition: Niggli.cc:120
bool meets_criteria_6(double compare_tol) const
|eta|<=A OR (eta==A, zeta<=2*ksi) OR (eta==-A, zeta==0)
Definition: Niggli.cc:156
bool meets_criteria_1(double compare_tol) const
A<=B OR (A==B, |ksi| <= |eta|)
Definition: Niggli.cc:73
bool meets_criteria_8(double compare_tol) const
Definition: Niggli.cc:204
double C() const
Square of lattice length c.
Definition: Niggli.cc:61
static std::vector< Eigen::Matrix3d > const & cell_invariant_transforms()
Definition: Niggli.cc:45
double A() const
Square of lattice length a.
Definition: Niggli.cc:57
bool is_niggli_type1(double compare_tol) const
True if all conditions except 4 are true.
Definition: Niggli.cc:224
bool is_niggli_type2(double compare_tol) const
True if all conditions except 3 are true.
Definition: Niggli.cc:228
void debug_criteria(double compare_tol) const
Definition: Niggli.cc:244
bool meets_criteria_2(double compare_tol) const
B<=C OR (B==C, |eta| <= |zeta|)
Definition: Niggli.cc:90
const Eigen::Matrix3d m_metrical_matrix
Transpose of initialization lattice dotted with itself.
Definition: Niggli.hh:123
const double m_scale_factor
Scaling factor for niggli comparisons.
Definition: Niggli.hh:126
double B() const
Square of lattice length b.
Definition: Niggli.cc:59
bool operator()(const Eigen::Matrix3d &LHS, const Eigen::Matrix3d &RHS) const
Definition: Niggli.hh:166
Index niggli_index(const Eigen::Matrix3d &test_lat_mat, double compare_tol)
Number of niggli criteria met.
Definition: Niggli.cc:267
bool standard_orientation_spatial_compare(const Eigen::Matrix3d &low_score_lat_mat, Eigen::Matrix3d &high_score_lat_mat, double compare_tol)
Compare the spatial orientation (ignoring matrix symmetry) and determine which one is oriented more s...
bool is_niggli(const Eigen::Matrix3d &test_lat_mat, double compare_tol)
Definition: Niggli.cc:271
Eigen::VectorXd spatial_unroll(const Eigen::Matrix3d &lat_mat, double compare_tol)
Generate a vector whose lexicographical value determines how well it's oriented in space.
Definition: Niggli.cc:354
Lattice niggli(const Lattice &in_lat, double compare_tol, bool keep_handedness=false)
Definition: Niggli.cc:330
bool standard_orientation_compare(const Eigen::Matrix3d &low_score_lat_mat, const Eigen::Matrix3d &high_score_lat_mat, double compare_tol)
Determine whether high_score has a more standard format than low_score.
Definition: Niggli.cc:427
Main CASM namespace.
Definition: APICommand.hh:8
Eigen::Matrix3d Matrix3d
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39
Eigen::VectorXd VectorXd