CASM  1.1.0
A Clusters Approach to Statistical Mechanics
Lattice.hh
Go to the documentation of this file.
1 #ifndef LATTICE_HH
2 #define LATTICE_HH
3 
4 #include <cmath>
5 #include <iostream>
6 #include <type_traits>
7 #include <vector>
8 
10 #include "casm/global/eigen.hh"
11 #include "casm/misc/Comparisons.hh"
12 
13 namespace CASM {
14 
15 namespace xtal {
16 
29 class Lattice : public Comparisons<CRTPBase<Lattice>> {
30  public:
31  typedef Eigen::Matrix3d::ColXpr LatVec;
32  typedef Eigen::Matrix3d::ConstColXpr ConstLatVec;
33 
34  Lattice(Eigen::Ref<const Eigen::Vector3d> const &vec1,
35  Eigen::Ref<const Eigen::Vector3d> const &vec2,
36  Eigen::Ref<const Eigen::Vector3d> const &vec3, double xtal_tol = TOL,
37  bool force = false);
38 
42  Lattice(Eigen::Ref<const Eigen::Matrix3d> const &lat_mat =
44  double xtal_tol = TOL, bool force = false);
45 
47  static Lattice fcc(double tol = TOL);
48 
50  static Lattice bcc(double tol = TOL);
51 
53  static Lattice cubic(double tol = TOL);
54 
56  static Lattice hexagonal(double tol = TOL);
57 
58  static std::vector<Eigen::Matrix3d> const &skew_transforms();
59 
61  LatVec operator[](Index i) { return m_lat_mat.col(i); }
62 
64  ConstLatVec operator[](Index i) const { return m_lat_mat.col(i); }
65 
66  std::tuple<LatVec, LatVec, LatVec> vectors() {
67  return std::make_tuple(m_lat_mat.col(0), m_lat_mat.col(1),
68  m_lat_mat.col(2));
69  }
70 
71  std::tuple<ConstLatVec, ConstLatVec, ConstLatVec> vectors() const {
72  return std::make_tuple(m_lat_mat.col(0), m_lat_mat.col(1),
73  m_lat_mat.col(2));
74  }
75 
78  Lattice scaled_lattice(double scale) const;
79 
82  double angle(Index i) const;
83 
85  double length(Index i) const;
86 
88  double volume() const { return lat_column_mat().determinant(); }
89 
95  Eigen::MatrixXd const &voronoi_table() const {
96  if (!m_voronoi_table.size()) {
98  }
99  return m_voronoi_table;
100  }
101 
104  double inner_voronoi_radius() const {
105  voronoi_table();
106  return m_inner_voronoi_radius;
107  }
108 
110  const Eigen::Matrix3d &lat_column_mat() const { return m_lat_mat; }
111 
118 
121  std::vector<int> calc_kpoints(std::vector<int> prim_kpoints,
122  Lattice prim_lat);
123 
125  Lattice reciprocal() const;
126 
128  double boxiness() const;
129 
133  Lattice reduced_cell() const;
134 
138  Lattice reduced_cell2() const;
139 
140  void print_voronoi_table(std::ostream &stream) const;
141 
143  double min_voronoi_radius() const;
144 
152  double max_voronoi_measure(const Eigen::Vector3d &pos,
153  Eigen::Vector3d &lattice_trans) const;
154 
161  int voronoi_number(const Eigen::Vector3d &pos) const;
162 
166  Eigen::Vector3i enclose_sphere(double radius) const;
167 
168  // TODO: Extract
169  void read(std::istream &stream);
170  // TODO: Extract
171  void print(std::ostream &stream, int _prec = 8) const;
172 
174  bool operator<(const Lattice &RHS) const;
175 
177  Lattice box(const Lattice &prim, const Lattice &scel,
178  bool verbose = false) const;
179 
182 
184  bool is_right_handed() const;
185 
188  Eigen::Vector3i millers(Eigen::Vector3d plane_normal) const;
189 
192  Lattice lattice_in_plane(Eigen::Vector3i millers,
193  int max_vol = 0) const; // John G 121030
194 
195  double tol() const { return m_tol; }
196 
197  void set_tol(double _tol) { m_tol = _tol; }
198 
199  private:
201 
203  bool _eq(const Lattice &RHS) const;
204 
206  void _generate_voronoi_table() const;
207 
208  mutable double m_inner_voronoi_radius;
210 
211  // Coordinate Conversion Matrices
212  // 0 is fractional to cartesion; 1 is cartesian to fractional
213  // use FRAC and CART globals to index
214  // e.g., coord[CART]=coord_trans[FRAC]*coord[FRAC]; //converts frac to cart
215  // Word to the wise: coord_trans[FRAC] is the matrix with columns equal to the
216  // lattice vectors
218 
219  double m_tol;
220 };
221 
229 double volume(const Lattice &lat);
230 
233 std::pair<bool, Eigen::Matrix3d> is_superlattice(const Lattice &scel,
234  const Lattice &unit,
235  double tol);
236 
237 std::istream &operator>>(std::istream &in, const Lattice &lattice_in);
238 
241 Lattice make_superduperlattice(const Lattice &lat1, const Lattice &lat2);
242 
248 Lattice replace_vector(const Lattice &lat, const Eigen::Vector3d &new_vector,
249  double tol);
250 
259  const Eigen::Ref<const Eigen::Matrix3d> &cart_mat, const Lattice &lat) {
260  return lat.inv_lat_column_mat() * cart_mat * lat.lat_column_mat();
261 }
262 
271  const Eigen::Ref<const Eigen::Matrix3d> &frac_mat, const Lattice &lat) {
272  return lat.lat_column_mat() * frac_mat * lat.inv_lat_column_mat();
273 }
274 
281 inline double volume(const Lattice &lat) {
282  return lat.lat_column_mat().determinant();
283 }
284 
286 template <typename IntegralType, int Options = 0>
288  const Lattice &lat,
289  const Eigen::Matrix<IntegralType, 3, 3, Options> &transf_mat) {
290  // TODO: Do this in more places that involve transformation matrices
291  static_assert(std::is_integral<IntegralType>::value,
292  "Transfomration matrix must be integer matrix");
293  return Lattice(Eigen::Matrix3d(lat.lat_column_mat()) *
294  transf_mat.template cast<double>(),
295  lat.tol());
296 }
297 
304  const Lattice &superlattice,
305  double tol);
306 
307 } // namespace xtal
308 } // namespace CASM
309 #endif
bool is_right_handed() const
Check if the lattice is right handed.
Definition: Lattice.cc:751
Eigen::Matrix3d::ConstColXpr ConstLatVec
Definition: Lattice.hh:32
Eigen::Matrix3d m_inv_lat_mat
Definition: Lattice.hh:217
void print(std::ostream &stream, int _prec=8) const
Definition: Lattice.cc:147
int voronoi_number(const Eigen::Vector3d &pos) const
Definition: Lattice.cc:362
Eigen::Vector3i millers(Eigen::Vector3d plane_normal) const
Definition: Lattice.cc:476
Lattice box(const Lattice &prim, const Lattice &scel, bool verbose=false) const
Return a lattice with diagonal matrix that fits around starting lattice.
std::vector< int > calc_kpoints(std::vector< int > prim_kpoints, Lattice prim_lat)
Definition: Lattice.cc:171
double length(Index i) const
Return length of i'th lattice vector.
Definition: Lattice.cc:127
const Eigen::Matrix3d & lat_column_mat() const
3x3 matrix with lattice vectors as its columne
Definition: Lattice.hh:110
ConstLatVec operator[](Index i) const
Get i'th lattice vector as column expression.
Definition: Lattice.hh:64
double m_inner_voronoi_radius
Definition: Lattice.hh:208
Eigen::Matrix3d::ColXpr LatVec
Definition: Lattice.hh:31
void set_tol(double _tol)
Definition: Lattice.hh:197
Eigen::Vector3i enclose_sphere(double radius) const
Definition: Lattice.cc:444
double inner_voronoi_radius() const
Radius of the largest sphere that fits wholly within the voronoi cell.
Definition: Lattice.hh:104
double max_voronoi_measure(const Eigen::Vector3d &pos, Eigen::Vector3d &lattice_trans) const
Definition: Lattice.cc:353
double tol() const
Definition: Lattice.hh:195
Lattice(Eigen::Ref< const Eigen::Vector3d > const &vec1, Eigen::Ref< const Eigen::Vector3d > const &vec2, Eigen::Ref< const Eigen::Vector3d > const &vec3, double xtal_tol=TOL, bool force=false)
Definition: Lattice.cc:12
std::tuple< LatVec, LatVec, LatVec > vectors()
Definition: Lattice.hh:66
static Lattice bcc(double tol=TOL)
Construct BCC primitive cell of unit volume.
Definition: Lattice.cc:97
Lattice lattice_in_plane(Eigen::Vector3i millers, int max_vol=0) const
Definition: Lattice.cc:496
bool _eq(const Lattice &RHS) const
Are lattice vectors identical for two lattices, within TOL.
Definition: Lattice.cc:747
std::tuple< ConstLatVec, ConstLatVec, ConstLatVec > vectors() const
Definition: Lattice.hh:71
double boxiness() const
Return boxiness factor directly proportional to volume/SA ratio.
Definition: Lattice.cc:167
Lattice reduced_cell2() const
Find the lattice vectors which give most compact unit cell Compactness is measured by how close lat_c...
Definition: Lattice.cc:229
Lattice reduced_cell() const
Find the lattice vectors which give most compact unit cell Compactness is measured by how close lat_c...
Definition: Lattice.cc:302
Eigen::MatrixXd m_voronoi_table
Definition: Lattice.hh:209
Lattice scaled_lattice(double scale) const
Return scaled copy of this lattice (Note: Volume will be scaled by scale^3)
Definition: Lattice.cc:123
double angle(Index i) const
Return angle between lattice vectors (*this)[(i+1)%3] and (*this)[(i+2)%3], in degrees.
Definition: Lattice.cc:132
static std::vector< Eigen::Matrix3d > const & skew_transforms()
Definition: Lattice.cc:68
double volume() const
Return signed volume of this lattice.
Definition: Lattice.hh:88
double min_voronoi_radius() const
Radius of largest sphere that totally fits inside the voronoi cell.
const Eigen::Matrix3d & inv_lat_column_mat() const
Inverse of Lattice::lat_column_mat() It is the transformation matrix 'C2F', such that f = C2F * c whe...
Definition: Lattice.hh:117
bool operator<(const Lattice &RHS) const
Compare two Lattice.
Definition: Lattice.cc:736
static Lattice cubic(double tol=TOL)
Construct simple cubic primitive cell of unit volume.
Definition: Lattice.cc:108
static Lattice hexagonal(double tol=TOL)
Construct cubic primitive cell of unit volume.
Definition: Lattice.cc:112
Eigen::Matrix3d m_lat_mat
Definition: Lattice.hh:217
Eigen::MatrixXd const & voronoi_table() const
Return voronoi table, which specifies outward-pointing normals of Lattice Voronoi cell....
Definition: Lattice.hh:95
Lattice & make_right_handed()
Flip c vector if it's on the wrong side of a-b plane – return (*this)
Definition: Lattice.cc:466
Lattice reciprocal() const
Return reciprocal lattice.
Definition: Lattice.cc:163
void print_voronoi_table(std::ostream &stream) const
LatVec operator[](Index i)
Get i'th lattice vector as column expression.
Definition: Lattice.hh:61
static Lattice fcc(double tol=TOL)
Construct FCC primitive cell of unit volume.
Definition: Lattice.cc:86
void _generate_voronoi_table() const
populate voronoi information.
Definition: Lattice.cc:382
void read(std::istream &stream)
Definition: Lattice.cc:136
Eigen::Matrix3d frac2cart(const Eigen::Ref< const Eigen::Matrix3d > &frac_mat, const Lattice &lat)
Returns 'cart_mat' which is transformation of 'frac_mat' if cart_vec_after = cart_mat*cart_vec then f...
Definition: Lattice.hh:270
double volume(const Lattice &lat)
Returns the volume of a Lattice.
Definition: Lattice.hh:281
std::istream & operator>>(std::istream &in, const Lattice &lattice_in)
Lattice make_superlattice(const Lattice &lat, const Eigen::Matrix< IntegralType, 3, 3, Options > &transf_mat)
Returns a super Lattice. Transformation matrix must be integer.
Definition: Lattice.hh:287
std::pair< bool, Eigen::Matrix3d > is_superlattice(const Lattice &scel, const Lattice &unit, double tol)
Check if scel is a superlattice of unitcell unit and some integer transformation matrix T.
Definition: Lattice.cc:836
Lattice make_superduperlattice(const Lattice &lat1, const Lattice &lat2)
returns Lattice that is smallest possible superlattice of both input Lattice
Definition: Lattice.cc:783
Eigen::Matrix3d cart2frac(const Eigen::Ref< const Eigen::Matrix3d > &cart_mat, const Lattice &lat)
Returns 'frac_mat' which is transformation of 'cart_mat' if cart_vec_after = cart_mat*cart_vec then f...
Definition: Lattice.hh:258
Eigen::Matrix3l transf_mat(const Lattice &prim_lat, const Lattice &super_lat)
IdentitySymRepBuilder Identity()
Eigen::Matrix3l make_transformation_matrix_to_super(const Lattice &tiling_unit, const Lattice &superlattice, double tol)
Definition: Lattice.cc:848
Lattice replace_vector(const Lattice &lat, const Eigen::Vector3d &new_vector, double tol)
Definition: Lattice.cc:812
Main CASM namespace.
Definition: APICommand.hh:8
Eigen::MatrixXd MatrixXd
Eigen::Matrix3d Matrix3d
const double TOL
Definition: definitions.hh:30
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39
Matrix< long int, 3, 3 > Matrix3l
Definition: eigen.hh:12
Implements other comparisons in terms of '<'.
Definition: Comparisons.hh:25