31 typedef Eigen::Matrix3d::ColXpr
LatVec;
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,
42 Lattice(Eigen::Ref<const Eigen::Matrix3d>
const &lat_mat =
44 double xtal_tol =
TOL,
bool force =
false);
66 std::tuple<LatVec, LatVec, LatVec>
vectors() {
71 std::tuple<ConstLatVec, ConstLatVec, ConstLatVec>
vectors()
const {
121 std::vector<int>
calc_kpoints(std::vector<int> prim_kpoints,
153 Eigen::Vector3d &lattice_trans)
const;
169 void read(std::istream &stream);
171 void print(std::ostream &stream,
int _prec = 8)
const;
178 bool verbose =
false)
const;
188 Eigen::Vector3i
millers(Eigen::Vector3d plane_normal)
const;
193 int max_vol = 0)
const;
259 const Eigen::Ref<const Eigen::Matrix3d> &cart_mat,
const Lattice &lat) {
271 const Eigen::Ref<const Eigen::Matrix3d> &frac_mat,
const Lattice &lat) {
286 template <
typename IntegralType,
int Options = 0>
289 const Eigen::Matrix<IntegralType, 3, 3, Options> &
transf_mat) {
291 static_assert(std::is_integral<IntegralType>::value,
292 "Transfomration matrix must be integer matrix");
bool is_right_handed() const
Check if the lattice is right handed.
Eigen::Matrix3d::ConstColXpr ConstLatVec
Eigen::Matrix3d m_inv_lat_mat
void print(std::ostream &stream, int _prec=8) const
int voronoi_number(const Eigen::Vector3d &pos) const
Eigen::Vector3i millers(Eigen::Vector3d plane_normal) const
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)
double length(Index i) const
Return length of i'th lattice vector.
const Eigen::Matrix3d & lat_column_mat() const
3x3 matrix with lattice vectors as its columne
ConstLatVec operator[](Index i) const
Get i'th lattice vector as column expression.
double m_inner_voronoi_radius
Eigen::Matrix3d::ColXpr LatVec
void set_tol(double _tol)
Eigen::Vector3i enclose_sphere(double radius) const
double inner_voronoi_radius() const
Radius of the largest sphere that fits wholly within the voronoi cell.
double max_voronoi_measure(const Eigen::Vector3d &pos, Eigen::Vector3d &lattice_trans) const
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)
std::tuple< LatVec, LatVec, LatVec > vectors()
static Lattice bcc(double tol=TOL)
Construct BCC primitive cell of unit volume.
Lattice lattice_in_plane(Eigen::Vector3i millers, int max_vol=0) const
bool _eq(const Lattice &RHS) const
Are lattice vectors identical for two lattices, within TOL.
std::tuple< ConstLatVec, ConstLatVec, ConstLatVec > vectors() const
double boxiness() const
Return boxiness factor directly proportional to volume/SA ratio.
Lattice reduced_cell2() const
Find the lattice vectors which give most compact unit cell Compactness is measured by how close lat_c...
Lattice reduced_cell() const
Find the lattice vectors which give most compact unit cell Compactness is measured by how close lat_c...
Eigen::MatrixXd m_voronoi_table
Lattice scaled_lattice(double scale) const
Return scaled copy of this lattice (Note: Volume will be scaled by scale^3)
double angle(Index i) const
Return angle between lattice vectors (*this)[(i+1)%3] and (*this)[(i+2)%3], in degrees.
static std::vector< Eigen::Matrix3d > const & skew_transforms()
double volume() const
Return signed volume of this lattice.
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...
bool operator<(const Lattice &RHS) const
Compare two Lattice.
static Lattice cubic(double tol=TOL)
Construct simple cubic primitive cell of unit volume.
static Lattice hexagonal(double tol=TOL)
Construct cubic primitive cell of unit volume.
Eigen::Matrix3d m_lat_mat
Eigen::MatrixXd const & voronoi_table() const
Return voronoi table, which specifies outward-pointing normals of Lattice Voronoi cell....
Lattice & make_right_handed()
Flip c vector if it's on the wrong side of a-b plane – return (*this)
Lattice reciprocal() const
Return reciprocal lattice.
void print_voronoi_table(std::ostream &stream) const
LatVec operator[](Index i)
Get i'th lattice vector as column expression.
static Lattice fcc(double tol=TOL)
Construct FCC primitive cell of unit volume.
void _generate_voronoi_table() const
populate voronoi information.
void read(std::istream &stream)
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...
double volume(const Lattice &lat)
Returns the volume of a Lattice.
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.
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.
Lattice make_superduperlattice(const Lattice &lat1, const Lattice &lat2)
returns Lattice that is smallest possible superlattice of both input Lattice
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...
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)
Lattice replace_vector(const Lattice &lat, const Eigen::Vector3d &new_vector, double tol)
INDEX_TYPE Index
For long integer indexing:
Matrix< long int, 3, 3 > Matrix3l
Implements other comparisons in terms of '<'.