33 typedef Eigen::Matrix3d::ColXpr
LatVec;
36 Lattice(
const Eigen::Vector3d &vec1,
const Eigen::Vector3d &vec2,
37 const Eigen::Vector3d &vec3);
41 Lattice(
const Eigen::Ref<const Eigen::Matrix3d> &lat_mat = Eigen::Matrix3d::Identity());
65 std::tuple<LatVec, LatVec, LatVec>
vectors() {
69 std::tuple<ConstLatVec, ConstLatVec, ConstLatVec>
vectors()
const {
131 template<
typename SymOpIterator,
typename SymOpOutputIterator>
132 SymOpOutputIterator
find_invariant_subgroup(SymOpIterator begin, SymOpIterator end, SymOpOutputIterator result,
double pg_tol =
TOL)
const;
172 template <
typename T>
189 double max_voronoi_measure(
const Eigen::Vector3d &pos, Eigen::Vector3d &lattice_trans)
const;
203 template<
typename CoordType,
typename CoordType2>
206 void read(std::istream &stream);
207 void print(std::ostream &stream,
int _prec = 8)
const;
237 Eigen::Vector3i
get_millers(Eigen::Vector3d plane_normal,
double tolerance =
TOL)
const;
243 void pg_converge(
double small_tol,
double large_tol,
double increment);
307 template<
typename Object,
typename OpIterator>
321 template<
typename LatIterator,
typename SymOpIterator>
324 SymOpIterator op_begin = SymOpIterator(),
325 SymOpIterator op_end = SymOpIterator());
361 template <
typename T>
void _generate_voronoi_table() const
populate voronoi information.
Eigen::Matrix3d m_inv_lat_mat
bool is_canonical(double tol=TOL) const
Check if Lattice is in the canonical form.
Eigen::Matrix3d::ConstColXpr ConstLatVec
Lattice superdupercell(const Lattice &lat1, const Lattice &lat2)
returns Lattice that is smallest possible supercell of both input Lattice
ConstLatVec operator[](Index i) const
Get i'th lattice vector as column expression.
void from_json(ClexDescription &desc, const jsonParser &json)
Implements other comparisons in terms of '<'.
void generate_point_group(SymGroup &point_group, double pg_tol=TOL) const
Populate.
std::tuple< ConstLatVec, ConstLatVec, ConstLatVec > vectors() const
Array< double > pg_converge(double large_tol)
double length(Index i) const
Return length of i'th lattice vector.
static Lattice hexagonal()
Construct cubic primitive cell of unit volume.
bool _eq(const Lattice &RHS) const
Are lattice vectors identical for two lattices, within TOL.
bool is_equivalent(const Lattice &RHS, double tol) const
Are two lattices the same, even if they have different lattice vectors.
jsonParser & to_json(const ClexDescription &desc, jsonParser &json)
Data structure for holding supercell enumeration properties.
Object copy_apply(const Transform &f, Object obj, Args &&...args)
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::MatrixXd const & voronoi_table() const
Return voronoi table, which specifies outward-pointing normals of Lattice Voronoi cell...
Lattice make_supercell(const Lattice &lat, const Eigen::Matrix3i &transf_mat)
Returns a super Lattice.
Lattice get_reduced_cell() const
Find the lattice vectors which give most compact unit cell Compactness is measured by how close lat_c...
void print(std::ostream &stream, int _prec=8) const
static Lattice bcc()
Construct BCC primitive cell of unit volume.
const Eigen::Matrix3d & lat_column_mat() const
3x3 matrix with lattice vectors as its columne
Eigen::Matrix3d::ColXpr LatVec
Lattice make_supercell(const Eigen::Matrix< T, 3, 3 > &trans_mat) const
make a supercell of this lattice. Equivalent to Lattice(lat_column_mat()*trans_mat) ...
Lattice get_lattice_in_plane(Eigen::Vector3i millers, int max_vol=20) const
Generates a lattice with vectors a and b parallel to the plane described by the miller indeces...
double m_inner_voronoi_radius
SymGroup is a collection of symmetry operations that satisfy the group property The symmetry operatio...
Eigen::MatrixXd m_voronoi_table
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...
Lattice(const Eigen::Vector3d &vec1, const Eigen::Vector3d &vec2, const Eigen::Vector3d &vec3)
static Lattice cubic()
Construct simple cubic primitive cell of unit volume.
SymOp is the Coordinate representation of a symmetry operation it keeps fraction (FRAC) and Cartesian...
LatVec operator[](Index i)
Get i'th lattice vector as column expression.
Array< CoordType > gridstruc_build(double max_radius, double min_radius, Array< CoordType > basis, CoordType2 lat_point)
Make a grid of lattice sites such that min_radius <= distance <= max_radius from. ...
void read(std::istream &stream)
EigenIndex Index
For long integer indexing:
SymOp to_canonical(double tol=TOL) const
Returns the operation that applied to *this returns the canonical form.
Lattice canonical_form(double tol=TOL) const
Returns the canonical equivalent Lattice, using the point group of the Lattice.
bool is_supercell_of(const Lattice &tile, Eigen::Matrix3d &multimat, double _tol=TOL) const
Matrix that relates two lattices (e.g., strain or slat)
double min_voronoi_radius() const
Radius of largest sphere that totally fits inside the voronoi cell.
Eigen::Vector3i get_millers(Eigen::Vector3d plane_normal, double tolerance=TOL) const
Given a normal vector, a Vector3 containing the miller indeces for the lattice is generated...
void print_voronoi_table(std::ostream &stream) const
double vol() const
Return signed volume of this lattice.
std::tuple< LatVec, LatVec, LatVec > vectors()
Lattice box(const Lattice &prim, const Lattice &scel, bool verbose=false) const
Return a lattice with diagonal matrix that fits around starting lattice.
SymOpOutputIterator find_invariant_subgroup(SymOpIterator begin, SymOpIterator end, SymOpOutputIterator result, double pg_tol=TOL) const
Output the SymOp that leave this lattice invariant.
std::istream & operator>>(std::istream &_in, std::vector< T > &vec)
Lattice scaled_lattice(double scale) const
Return scaled copy of this lattice (Note: Volume will be scaled by scale^3)
bool operator<(const Lattice &RHS) const
Compare two Lattice.
Lattice get_reciprocal() const
Return reciprocal lattice.
bool is_right_handed() const
Check if the lattice is right handed.
double inner_voronoi_radius() const
Radius of the largest sphere that fits wholly within the voronoi cell.
Eigen::Vector3i enclose_sphere(double radius) const
double max_voronoi_measure(const Eigen::Vector3d &pos, Eigen::Vector3d &lattice_trans) const
void symmetrize(const SymGroup &relaxed_pg)
Force this lattice to have symmetry of group.
Lattice & make_right_handed()
Flip c vector if it's on the wrong side of a-b plane – return (*this)
std::pair< bool, Eigen::MatrixXi > is_supercell(const Lattice &scel, const Lattice &unit, double tol)
Check if scel is a supercell of unitcell unit and some integer transformation matrix T...
double volume(const Lattice &lat)
Returns the volume of a Lattice.
Array< int > calc_kpoints(Array< int > prim_kpoints, Lattice prim_lat)
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 angle(Index i) const
Return angle between lattice vectors (*this)[(i+1)%3] and (*this)[(i+2)%3], in degrees.
Lattice replace_vector(const Lattice &lat, const Eigen::Vector3d &new_vector, double tol)
SymOp from_canonical(double tol=TOL) const
Returns the operation that applied to the canonical form returns *this.
static Lattice fcc()
Construct FCC primitive cell of unit volume.
Object & apply(const Transform &f, Object &obj, Args &&...args)
Eigen::Matrix3d m_lat_mat
int voronoi_number(const Eigen::Vector3d &pos) const
void generate_supercells(Array< Lattice > &supercell, const SymGroup &effective_pg, const ScelEnumProps &enum_props) const
Populate.