CASM  1.1.0
A Clusters Approach to Statistical Mechanics
CASM::xtal::Lattice Class Reference

#include <Lattice.hh>

+ Inheritance diagram for CASM::xtal::Lattice:

Detailed Description

Definition at line 29 of file Lattice.hh.

Public Types

typedef Eigen::Matrix3d::ColXpr LatVec
 
typedef Eigen::Matrix3d::ConstColXpr ConstLatVec
 
typedef Base::MostDerived MostDerived
 

Public Member Functions

 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)
 
 Lattice (Eigen::Ref< const Eigen::Matrix3d > const &lat_mat=Eigen::Matrix3d::Identity(), double xtal_tol=TOL, bool force=false)
 
LatVec operator[] (Index i)
 Get i'th lattice vector as column expression. More...
 
ConstLatVec operator[] (Index i) const
 Get i'th lattice vector as column expression. More...
 
std::tuple< LatVec, LatVec, LatVecvectors ()
 
std::tuple< ConstLatVec, ConstLatVec, ConstLatVecvectors () const
 
Lattice scaled_lattice (double scale) const
 Return scaled copy of this lattice (Note: Volume will be scaled by scale^3) More...
 
double angle (Index i) const
 Return angle between lattice vectors (*this)[(i+1)%3] and (*this)[(i+2)%3], in degrees. More...
 
double length (Index i) const
 Return length of i'th lattice vector. More...
 
double volume () const
 Return signed volume of this lattice. More...
 
Eigen::MatrixXd const & voronoi_table () const
 Return voronoi table, which specifies outward-pointing normals of Lattice Voronoi cell. outward-pointing normals are given as rows of the matrix, and defined such that if (voronoi_table()*coord.cart()).maxCoeff()>1, then 'coord' is outside of the voronoi cell. More...
 
double inner_voronoi_radius () const
 Radius of the largest sphere that fits wholly within the voronoi cell. More...
 
const Eigen::Matrix3d & lat_column_mat () const
 3x3 matrix with lattice vectors as its columne More...
 
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 where 'c' is Cartesian coordinate/vector and 'f' and 'f' is fractional coordinate/vector. More...
 
std::vector< int > calc_kpoints (std::vector< int > prim_kpoints, Lattice prim_lat)
 
Lattice reciprocal () const
 Return reciprocal lattice. More...
 
double boxiness () const
 Return boxiness factor directly proportional to volume/SA ratio. More...
 
Lattice reduced_cell () const
 Find the lattice vectors which give most compact unit cell Compactness is measured by how close lat_column_mat().transpose()*lat_column_mat() is to a diagonal matrix. More...
 
Lattice reduced_cell2 () const
 Find the lattice vectors which give most compact unit cell Compactness is measured by how close lat_column_mat().transpose()*lat_column_mat() is to a diagonal matrix. More...
 
void print_voronoi_table (std::ostream &stream) const
 
double min_voronoi_radius () const
 Radius of largest sphere that totally fits inside the voronoi cell. More...
 
double max_voronoi_measure (const Eigen::Vector3d &pos, Eigen::Vector3d &lattice_trans) const
 
int voronoi_number (const Eigen::Vector3d &pos) const
 
Eigen::Vector3i enclose_sphere (double radius) const
 
void read (std::istream &stream)
 
void print (std::ostream &stream, int _prec=8) const
 
bool operator< (const Lattice &RHS) const
 Compare two Lattice. More...
 
Lattice box (const Lattice &prim, const Lattice &scel, bool verbose=false) const
 Return a lattice with diagonal matrix that fits around starting lattice. More...
 
Latticemake_right_handed ()
 Flip c vector if it's on the wrong side of a-b plane – return (*this) More...
 
bool is_right_handed () const
 Check if the lattice is right handed. More...
 
Eigen::Vector3i millers (Eigen::Vector3d plane_normal) const
 
Lattice lattice_in_plane (Eigen::Vector3i millers, int max_vol=0) const
 
double tol () const
 
void set_tol (double _tol)
 
bool operator> (const MostDerived &B) const
 
bool operator<= (const MostDerived &B) const
 
bool operator>= (const MostDerived &B) const
 
bool operator== (const MostDerived &B) const
 
bool operator!= (const MostDerived &B) const
 

Static Public Member Functions

static Lattice fcc (double tol=TOL)
 Construct FCC primitive cell of unit volume. More...
 
static Lattice bcc (double tol=TOL)
 Construct BCC primitive cell of unit volume. More...
 
static Lattice cubic (double tol=TOL)
 Construct simple cubic primitive cell of unit volume. More...
 
static Lattice hexagonal (double tol=TOL)
 Construct cubic primitive cell of unit volume. More...
 
static std::vector< Eigen::Matrix3d > const & skew_transforms ()
 

Protected Member Functions

bool eq_impl (const MostDerived &B) const
 
bool ne_impl (const MostDerived &B) const
 

Private Member Functions

bool _eq (const Lattice &RHS) const
 Are lattice vectors identical for two lattices, within TOL. More...
 
void _generate_voronoi_table () const
 populate voronoi information. More...
 

Private Attributes

friend Comparisons< Lattice >
 
double m_inner_voronoi_radius
 
Eigen::MatrixXd m_voronoi_table
 
Eigen::Matrix3d m_lat_mat
 
Eigen::Matrix3d m_inv_lat_mat
 
double m_tol
 

Related Functions

(Note that these are not member functions.)

Lattice replace_vector (const Lattice &lat, const Eigen::Vector3d &new_vector, double tol)
 Returns a minimum volume Lattice obtainable by replacing one Lattice vector. More...
 
Eigen::Matrix3i canonical_hnf (const Eigen::Matrix3i &T, const SymOpVector &effective_pg, const Lattice &ref_lattice)
 Return canonical hermite normal form of the supercell matrix. More...
 

Member Typedef Documentation

◆ ConstLatVec

typedef Eigen::Matrix3d::ConstColXpr CASM::xtal::Lattice::ConstLatVec

Definition at line 32 of file Lattice.hh.

◆ LatVec

typedef Eigen::Matrix3d::ColXpr CASM::xtal::Lattice::LatVec

Definition at line 31 of file Lattice.hh.

◆ MostDerived

template<typename Base >
typedef Base::MostDerived notstd::Comparisons< Base >::MostDerived
inherited

Definition at line 27 of file Comparisons.hh.

Constructor & Destructor Documentation

◆ Lattice() [1/2]

CASM::xtal::Lattice::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 at line 12 of file Lattice.cc.

◆ Lattice() [2/2]

CASM::xtal::Lattice::Lattice ( Eigen::Ref< const Eigen::Matrix3d > const &  lat_mat = Eigen::Matrix3d::Identity(),
double  xtal_tol = TOL,
bool  force = false 
)

Construct Lattice from a matrix of lattice vectors, where lattice vectors are columns (e.g., lat_mat is equivalent to coord_trans[FRAC])

Construct Lattice from a matrix of lattice vectors, where lattice vectors are columns (e.g., lat_mat is equivalent to lat_column_mat())

Definition at line 76 of file Lattice.cc.

Member Function Documentation

◆ _eq()

bool CASM::xtal::Lattice::_eq ( const Lattice RHS) const
private

Are lattice vectors identical for two lattices, within TOL.

Are lattice vectors identical for two lattices.

Definition at line 747 of file Lattice.cc.

◆ _generate_voronoi_table()

void CASM::xtal::Lattice::_generate_voronoi_table ( ) const
private

populate voronoi information.

Definition at line 382 of file Lattice.cc.

◆ angle()

double CASM::xtal::Lattice::angle ( Index  i) const

Return angle between lattice vectors (*this)[(i+1)%3] and (*this)[(i+2)%3], in degrees.

Definition at line 132 of file Lattice.cc.

◆ bcc()

Lattice CASM::xtal::Lattice::bcc ( double  tol = TOL)
static

Construct BCC primitive cell of unit volume.

Definition at line 97 of file Lattice.cc.

◆ box()

Lattice CASM::xtal::Lattice::box ( const Lattice prim,
const Lattice scel,
bool  verbose = false 
) const

Return a lattice with diagonal matrix that fits around starting lattice.

◆ boxiness()

double CASM::xtal::Lattice::boxiness ( ) const

Return boxiness factor directly proportional to volume/SA ratio.

Definition at line 167 of file Lattice.cc.

◆ calc_kpoints()

std::vector< int > CASM::xtal::Lattice::calc_kpoints ( std::vector< int >  prim_kpoints,
Lattice  prim_lat 
)

Calculates the kpoint mesh for a superlattice lattice given the kpoint mesh for the primitive lattice

Definition at line 171 of file Lattice.cc.

◆ cubic()

Lattice CASM::xtal::Lattice::cubic ( double  tol = TOL)
static

Construct simple cubic primitive cell of unit volume.

Definition at line 108 of file Lattice.cc.

◆ enclose_sphere()

Eigen::Vector3i CASM::xtal::Lattice::enclose_sphere ( double  radius) const

Gives a scaling of the lattice vectors such that after scaling, The eight parallelipipeds touching the origin enclose a sphere of given radius

This function, given a linearly independent set of lattice vectors, finds the dimensions along the unit cell vectors such that a sphere of given radius fits within a uniform grid of 2dim[1]x2dim[2]x2dim[3] lattice points centered at the origin.

The algorithm works by getting the normal (e.g. n1) to each pair of lattice vectors (e.g. a2, a3), scaling this normal to have length radius and then projecting this normal parallel to the a2,a3 plane onto the remaining lattice vector a1. This will tell us the number of a1 vectors needed to make a grid to encompass the sphere.

Definition at line 444 of file Lattice.cc.

◆ eq_impl()

template<typename Base >
bool notstd::Comparisons< Base >::eq_impl ( const MostDerived B) const
inlineprotectedinherited

Definition at line 40 of file Comparisons.hh.

◆ fcc()

Lattice CASM::xtal::Lattice::fcc ( double  tol = TOL)
static

Construct FCC primitive cell of unit volume.

Definition at line 86 of file Lattice.cc.

◆ hexagonal()

Lattice CASM::xtal::Lattice::hexagonal ( double  tol = TOL)
static

Construct cubic primitive cell of unit volume.

Definition at line 112 of file Lattice.cc.

◆ inner_voronoi_radius()

double CASM::xtal::Lattice::inner_voronoi_radius ( ) const
inline

Radius of the largest sphere that fits wholly within the voronoi cell.

Definition at line 104 of file Lattice.hh.

◆ inv_lat_column_mat()

const Eigen::Matrix3d& CASM::xtal::Lattice::inv_lat_column_mat ( ) const
inline

Inverse of Lattice::lat_column_mat() It is the transformation matrix 'C2F', such that f = C2F * c where 'c' is Cartesian coordinate/vector and 'f' and 'f' is fractional coordinate/vector.

Definition at line 117 of file Lattice.hh.

◆ is_right_handed()

bool CASM::xtal::Lattice::is_right_handed ( ) const

Check if the lattice is right handed.

Definition at line 751 of file Lattice.cc.

◆ lat_column_mat()

const Eigen::Matrix3d& CASM::xtal::Lattice::lat_column_mat ( ) const
inline

3x3 matrix with lattice vectors as its columne

Definition at line 110 of file Lattice.hh.

◆ lattice_in_plane()

Lattice CASM::xtal::Lattice::lattice_in_plane ( Eigen::Vector3i  millers,
int  max_vol = 0 
) const

Generates a lattice with vectors a and b parallel to the plane described by the miller indeces

Using miller indeces, the intercept of the plane is calculated. All intercepts are multiplied by a factor such that their values are integers. This corresponds to the plane being shifted so that it lands on three lattice sites. Using these lattice sites new vectors A and B that lie on the desired plane can be constructed. The choice of vector C is somewhat arbitrary. lattice_in_plane will first construct C to be the shortest most orthogonal vector possible. The user is then given the option to make C more orthogonal to A and B in exchange for a greater length. The returned lattice will be the one containing the most orthogonal C that is still smaller than max_vol

Definition at line 496 of file Lattice.cc.

◆ length()

double CASM::xtal::Lattice::length ( Index  i) const

Return length of i'th lattice vector.

Definition at line 127 of file Lattice.cc.

◆ make_right_handed()

Lattice & CASM::xtal::Lattice::make_right_handed ( )

Flip c vector if it's on the wrong side of a-b plane – return (*this)

A lattice is considered right handed when the determinant of the lattice vector matrix is positive. This routine will flip the sign of all the lattice vectors if it finds that the determinant is negative

Definition at line 466 of file Lattice.cc.

◆ max_voronoi_measure()

double CASM::xtal::Lattice::max_voronoi_measure ( const Eigen::Vector3d &  pos,
Eigen::Vector3d &  lattice_trans 
) const

Given cartesian vector 'pos', returns double v_dist and populates 'lattice_trans'. v_dist is the fractional number of half lattice-planes along 'lattice_trans' between 'pos' and the lattice plane that passes through the origin. lattice_trans is a lattice translation that is normal to a face of the voronoi cell and translating pos to 'pos-lattice_trans' will yield a translationally-equivalent vector between the -1 and +1 half-plane. v_dist is maximized over all faces of the voronoi cell.

Definition at line 353 of file Lattice.cc.

◆ millers()

Eigen::Vector3i CASM::xtal::Lattice::millers ( Eigen::Vector3d  plane_normal) const

Given a normal vector, a Vector3 containing the miller indeces for the lattice is generated

Definition at line 476 of file Lattice.cc.

◆ min_voronoi_radius()

double CASM::xtal::Lattice::min_voronoi_radius ( ) const

Radius of largest sphere that totally fits inside the voronoi cell.

◆ ne_impl()

template<typename Base >
bool notstd::Comparisons< Base >::ne_impl ( const MostDerived B) const
inlineprotectedinherited

Definition at line 44 of file Comparisons.hh.

◆ operator!=()

template<typename Base >
bool notstd::Comparisons< Base >::operator!= ( const MostDerived B) const
inlineinherited

Definition at line 37 of file Comparisons.hh.

◆ operator<()

bool CASM::xtal::Lattice::operator< ( const Lattice RHS) const

Compare two Lattice.

  • First compares is_niggli(*this, TOL) with is_niggli(RHS, TOL)
  • Then compares via standard_orientation_compare(this->lat_column_mat(), RHS.lat_column_mat(), TOL)

Definition at line 736 of file Lattice.cc.

◆ operator<=()

template<typename Base >
bool notstd::Comparisons< Base >::operator<= ( const MostDerived B) const
inlineinherited

Definition at line 31 of file Comparisons.hh.

◆ operator==()

template<typename Base >
bool notstd::Comparisons< Base >::operator== ( const MostDerived B) const
inlineinherited

Definition at line 35 of file Comparisons.hh.

◆ operator>()

template<typename Base >
bool notstd::Comparisons< Base >::operator> ( const MostDerived B) const
inlineinherited

Definition at line 29 of file Comparisons.hh.

◆ operator>=()

template<typename Base >
bool notstd::Comparisons< Base >::operator>= ( const MostDerived B) const
inlineinherited

Definition at line 33 of file Comparisons.hh.

◆ operator[]() [1/2]

LatVec CASM::xtal::Lattice::operator[] ( Index  i)
inline

Get i'th lattice vector as column expression.

Definition at line 61 of file Lattice.hh.

◆ operator[]() [2/2]

ConstLatVec CASM::xtal::Lattice::operator[] ( Index  i) const
inline

Get i'th lattice vector as column expression.

Definition at line 64 of file Lattice.hh.

◆ print()

void CASM::xtal::Lattice::print ( std::ostream &  stream,
int  _prec = 8 
) const

Definition at line 147 of file Lattice.cc.

◆ print_voronoi_table()

void CASM::xtal::Lattice::print_voronoi_table ( std::ostream &  stream) const

◆ read()

void CASM::xtal::Lattice::read ( std::istream &  stream)

Definition at line 136 of file Lattice.cc.

◆ reciprocal()

Lattice CASM::xtal::Lattice::reciprocal ( ) const

Return reciprocal lattice.

Definition at line 163 of file Lattice.cc.

◆ reduced_cell()

Lattice CASM::xtal::Lattice::reduced_cell ( ) const

Find the lattice vectors which give most compact unit cell Compactness is measured by how close lat_column_mat().transpose()*lat_column_mat() is to a diagonal matrix.

This function finds the reduced cell from the given primitive cell.

This implementation is the LLL algorithm as laid out by Hoffstein, Jeffrey; Pipher, Jill; Silverman, J.H. (2008). An Introduction to Mathematical Cryptography. Springer. ISBN 978-0-387-77993-5. Some corrections have been made from Silverman, Joseph. "Introduction to Mathematical Cryptography Errata" (PDF). Brown University Mathematics Dept. Retrieved 5 May 2015.

Definition at line 302 of file Lattice.cc.

◆ reduced_cell2()

Lattice CASM::xtal::Lattice::reduced_cell2 ( ) const

Find the lattice vectors which give most compact unit cell Compactness is measured by how close lat_column_mat().transpose()*lat_column_mat() is to a diagonal matrix.

This function finds the reduced cell from the given primitive cell.

First, the translation vectors of the primitive cell are determined, then reduced. The Bravais lattice is determined using Niggli's transformations. The method used to determine the basis vectors are to find body and face diagonals across the cell. The vectors satisfying these main conditions can be found by ensuring: i) the shortest face diagonal is not shorter than the longest edge of the same face ii) the shortest body diagonal is not shorter than the longest edge of the cell

This is accomplished by enumerating the skew matrices (identity matrix with one of the off-diagonal terms =1 or =-1). This will replace one of the lattice vectors with the face diagonal corresponding to that linear combination. If the length of any of the new skewed lattice vectors is shorter than the length of the original cell, it is replaced. All skew matrices have determinants equal to 1, which preserves the volume of the original cell.

Definition at line 229 of file Lattice.cc.

◆ scaled_lattice()

Lattice CASM::xtal::Lattice::scaled_lattice ( double  scale) const

Return scaled copy of this lattice (Note: Volume will be scaled by scale^3)

Definition at line 123 of file Lattice.cc.

◆ set_tol()

void CASM::xtal::Lattice::set_tol ( double  _tol)
inline

Definition at line 197 of file Lattice.hh.

◆ skew_transforms()

std::vector< Eigen::Matrix3d > const & CASM::xtal::Lattice::skew_transforms ( )
static

Definition at line 68 of file Lattice.cc.

◆ tol()

double CASM::xtal::Lattice::tol ( ) const
inline

Definition at line 195 of file Lattice.hh.

◆ vectors() [1/2]

std::tuple<LatVec, LatVec, LatVec> CASM::xtal::Lattice::vectors ( )
inline

Definition at line 66 of file Lattice.hh.

◆ vectors() [2/2]

std::tuple<ConstLatVec, ConstLatVec, ConstLatVec> CASM::xtal::Lattice::vectors ( ) const
inline

Definition at line 71 of file Lattice.hh.

◆ volume()

double CASM::xtal::Lattice::volume ( ) const
inline

Return signed volume of this lattice.

Definition at line 88 of file Lattice.hh.

◆ voronoi_number()

int CASM::xtal::Lattice::voronoi_number ( const Eigen::Vector3d &  pos) const

return number of voronoi cell faces that 'pos' is on, within +/- TOL 0 indicates that 'pos' is within the voronoi cell, and the origin is the nearest lattice site -1 indicates that 'pos' is outside the voronoi cell, and there is a lattice site closer than the origin Values of 1<=n<=7 indicate that there n lattice sites equally as close as the origin

Definition at line 362 of file Lattice.cc.

◆ voronoi_table()

Eigen::MatrixXd const& CASM::xtal::Lattice::voronoi_table ( ) const
inline

Return voronoi table, which specifies outward-pointing normals of Lattice Voronoi cell. outward-pointing normals are given as rows of the matrix, and defined such that if (voronoi_table()*coord.cart()).maxCoeff()>1, then 'coord' is outside of the voronoi cell.

Definition at line 95 of file Lattice.hh.

Member Data Documentation

◆ Comparisons< Lattice >

friend CASM::xtal::Lattice::Comparisons< Lattice >
private

Definition at line 200 of file Lattice.hh.

◆ m_inner_voronoi_radius

double CASM::xtal::Lattice::m_inner_voronoi_radius
mutableprivate

Definition at line 208 of file Lattice.hh.

◆ m_inv_lat_mat

Eigen::Matrix3d CASM::xtal::Lattice::m_inv_lat_mat
private

Definition at line 217 of file Lattice.hh.

◆ m_lat_mat

Eigen::Matrix3d CASM::xtal::Lattice::m_lat_mat
private

Definition at line 217 of file Lattice.hh.

◆ m_tol

double CASM::xtal::Lattice::m_tol
private

Definition at line 219 of file Lattice.hh.

◆ m_voronoi_table

Eigen::MatrixXd CASM::xtal::Lattice::m_voronoi_table
mutableprivate

Definition at line 209 of file Lattice.hh.


The documentation for this class was generated from the following files: