CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
CASM::PrimNeighborList Class Reference

#include <NeighborList.hh>

Detailed Description

The PrimNeighborList gives the coordinates of UnitCell that are neighbors of the origin UnitCell.

  • The PrimNeighborList is constructed with a weighting matrix, M, that defines the shape of neighborhood
  • The canonical order of neighboring UnitCell is obtained by lexicographically sorting [r, i, j, k], where r = (i,j,k).transpose() * W * (i,j,k).
  • The canonical order of UnitCellCoord is obtained by lexicographically sorting [r, i, j, k, b], where r = (i,j,k).transpose() * W * (i,j,k).
  • The PrimNeighborList can be expanded as needed to increase the range of the neighborhood without affecting the order of the neighbors.

Definition at line 24 of file NeighborList.hh.

Public Types

typedef Eigen::Matrix3l Matrix3Type
 
typedef Eigen::Matrix< long,
Eigen::Dynamic, 1 > 
VectorXType
 
typedef std::set< UnitCell,
std::function< bool(UnitCell,
UnitCell)> > 
NeighborSet
 
typedef NeighborSet::size_type size_type
 
typedef NeighborSet::const_iterator const_iterator
 
typedef std::set< int > SublatIndices
 

Public Member Functions

Matrix3Type weight_matrix () const
 Return the weighting matrix W used to define the canonical order. More...
 
void expand (UnitCellCoord uccoord)
 Expand the neighbor list to include the given UnitCellCoord. More...
 
size_type size () const
 size of the neighborhood of unit cells More...
 
const_iterator begin () const
 const_iterator over the neighborhood of unit cells More...
 
const_iterator end () const
 const_iterator over the neighborhood of unit cells More...
 
const_iterator cbegin () const
 const_iterator over the neighborhood of unit cells More...
 
const_iterator cend () const
 const_iterator over the neighborhood of unit cells More...
 
const SublatIndicessublat_indices () const
 pair of const_iterators over a range of indices of sublattices to include More...
 
std::unique_ptr< PrimNeighborListclone () const
 Clone. More...
 
template<typename SublatIterator >
 PrimNeighborList (const Matrix3Type &W, SublatIterator begin, SublatIterator end)
 Constructor, specifying the weighting matrix to use and the indices of sublattices to include. More...
 
template<typename UnitCellCoordIterator >
void expand (UnitCellCoordIterator begin, UnitCellCoordIterator end)
 Expand the neighbor list to include the given range of UnitCellCoord. More...
 

Static Public Member Functions

static Matrix3Type make_weight_matrix (const Eigen::Matrix3d lat_column_mat, Index max_element_value, double tol)
 Returns a NeighborList weighting matrix appropriate for a particular lattice. More...
 

Public Attributes

EIGEN_MAKE_ALIGNED_OPERATOR_NEW
typedef Index 
Scalar
 

Private Member Functions

void _expand (Scalar prev_range)
 Ensure that all intermediate UnitCell are included in our neighborhood. More...
 
Scalar _dist (const UnitCell &A) const
 Calculate A.transpose()*M*A. More...
 
VectorXType _add_dist (const UnitCell &A) const
 Return [r, i, j, k], where r = _dist(A) More...
 
bool _compare_unitcell (const UnitCell &A, const UnitCell &B) const
 Convert [i,j,k] -> [r,i,j,k] and then lexicographically compare. More...
 

Static Private Member Functions

static bool _compare_vec (const VectorXType &A, const VectorXType &B)
 Lexicographical comparison. More...
 

Private Attributes

Matrix3Type m_W
 Weighting matrix. More...
 
Eigen::MatrixXd m_Uinv
 Cholesky decomposition: W = U.transpose()*U, useful for getting bounding box. More...
 
NeighborSet m_neighborhood
 the neighborhood of [i, j, k] More...
 
Scalar m_range
 the neighborhood, m_neighborhood, contains all UnitCell with r <= m_range More...
 
SublatIndices m_sublat_indices
 the indices of sublattices that should be included More...
 

Member Typedef Documentation

typedef NeighborSet::const_iterator CASM::PrimNeighborList::const_iterator

Definition at line 35 of file NeighborList.hh.

Definition at line 33 of file NeighborList.hh.

typedef NeighborSet::size_type CASM::PrimNeighborList::size_type

Definition at line 34 of file NeighborList.hh.

Definition at line 36 of file NeighborList.hh.

typedef Eigen::Matrix<long, Eigen::Dynamic, 1> CASM::PrimNeighborList::VectorXType

Definition at line 32 of file NeighborList.hh.

Constructor & Destructor Documentation

template<typename SublatIterator >
CASM::PrimNeighborList::PrimNeighborList ( const Matrix3Type W,
SublatIterator  begin,
SublatIterator  end 
)

Constructor, specifying the weighting matrix to use and the indices of sublattices to include.

Constructor.

Parameters
WA positive-definite weighting matrix for calculating the weighted squared L2 norm of the unitcell indices.
  • The canonical order of UnitCell is obtained by lexicographically sorting [r, i, j, k], where r = (i,j,k).transpose() * W * (i,j,k).
  • The canonical order of UnitCellCoord is obtained by lexicographically sorting [r, i, j, k, b], where r = (i,j,k).transpose() * W * (i,j,k).

Definition at line 188 of file NeighborList.hh.

Member Function Documentation

PrimNeighborList::VectorXType CASM::PrimNeighborList::_add_dist ( const UnitCell A) const
private

Return [r, i, j, k], where r = _dist(A)

Definition at line 112 of file NeighborList.cc.

bool CASM::PrimNeighborList::_compare_unitcell ( const UnitCell A,
const UnitCell B 
) const
private

Convert [i,j,k] -> [r,i,j,k] and then lexicographically compare.

Definition at line 120 of file NeighborList.cc.

bool CASM::PrimNeighborList::_compare_vec ( const VectorXType A,
const VectorXType B 
)
staticprivate

Lexicographical comparison.

Definition at line 125 of file NeighborList.cc.

PrimNeighborList::Scalar CASM::PrimNeighborList::_dist ( const UnitCell A) const
private

Calculate A.transpose()*M*A.

Definition at line 107 of file NeighborList.cc.

void CASM::PrimNeighborList::_expand ( Scalar  prev_range)
private

Ensure that all intermediate UnitCell are included in our neighborhood.

Definition at line 29 of file NeighborList.cc.

PrimNeighborList::const_iterator CASM::PrimNeighborList::begin ( ) const

const_iterator over the neighborhood of unit cells

Definition at line 82 of file NeighborList.cc.

PrimNeighborList::const_iterator CASM::PrimNeighborList::cbegin ( ) const

const_iterator over the neighborhood of unit cells

Definition at line 92 of file NeighborList.cc.

PrimNeighborList::const_iterator CASM::PrimNeighborList::cend ( ) const

const_iterator over the neighborhood of unit cells

Definition at line 97 of file NeighborList.cc.

std::unique_ptr< PrimNeighborList > CASM::PrimNeighborList::clone ( ) const

Clone.

Definition at line 171 of file NeighborList.cc.

PrimNeighborList::const_iterator CASM::PrimNeighborList::end ( ) const

const_iterator over the neighborhood of unit cells

Definition at line 87 of file NeighborList.cc.

void CASM::PrimNeighborList::expand ( UnitCellCoord  uccoord)

Expand the neighbor list to include the given UnitCellCoord.

Definition at line 14 of file NeighborList.cc.

template<typename UnitCellCoordIterator >
void CASM::PrimNeighborList::expand ( UnitCellCoordIterator  begin,
UnitCellCoordIterator  end 
)

Expand the neighbor list to include the given range of UnitCellCoord.

  • The canonical order of UnitCellCoord is obtained by lexicographically sorting [r, i, j, k, b], where r = (i,j,k).transpose() * W * (i,j,k).
  • The canonical order of UnitCell is obtained by lexicographically sorting [r, i, j, k]
  • After determining the sorted order of the requested UnitCell, all intermediate UnitCell will also be added so that indices into the neighbor list are always constant and in canonical order

Definition at line 221 of file NeighborList.hh.

PrimNeighborList::Matrix3Type CASM::PrimNeighborList::make_weight_matrix ( const Eigen::Matrix3d  lat_column_mat,
Index  max_element_value,
double  tol 
)
static

Returns a NeighborList weighting matrix appropriate for a particular lattice.

Returns the integer weight matrix, M, closest to lat.transpose()*lat with all elements less than max_element_value.

  • C = L*F, where L is the lattice vectors as a column matrix
  • Want r = F.transpose*M*F such that r = C.transpose*C, where F is fractional coordinates, and C is cartesian coordinates.
  • Find F.t*M*F = F.t*L.t*L*F -> W ~ L.tranpose*L
  • Approch, calculate L.tranpose*L, then divide by minimum element -> M'
  • consider integer n=1,2,...,n_last for all n such that all elements of n*M' are less than 'max_element_value'
  • return the first n*M' that is integer. If none are integer, return n_last*M'

Definition at line 142 of file NeighborList.cc.

PrimNeighborList::size_type CASM::PrimNeighborList::size ( ) const

size of the neighborhood of unit cells

Definition at line 77 of file NeighborList.cc.

const PrimNeighborList::SublatIndices & CASM::PrimNeighborList::sublat_indices ( ) const

pair of const_iterators over a range of indices of sublattices to include

const_iterator over the neighborhood of unit cells

Definition at line 102 of file NeighborList.cc.

PrimNeighborList::Matrix3Type CASM::PrimNeighborList::weight_matrix ( ) const

Return the weighting matrix W used to define the canonical order.

Return the weighting matrix used to define the canonical order.

Definition at line 9 of file NeighborList.cc.

Member Data Documentation

NeighborSet CASM::PrimNeighborList::m_neighborhood
private

the neighborhood of [i, j, k]

Definition at line 108 of file NeighborList.hh.

Scalar CASM::PrimNeighborList::m_range
private

the neighborhood, m_neighborhood, contains all UnitCell with r <= m_range

Definition at line 111 of file NeighborList.hh.

SublatIndices CASM::PrimNeighborList::m_sublat_indices
private

the indices of sublattices that should be included

Definition at line 114 of file NeighborList.hh.

Eigen::MatrixXd CASM::PrimNeighborList::m_Uinv
private

Cholesky decomposition: W = U.transpose()*U, useful for getting bounding box.

Definition at line 105 of file NeighborList.hh.

Matrix3Type CASM::PrimNeighborList::m_W
private

Weighting matrix.

  • The canonical order of UnitCellCoord is obtained by lexicographically sorting [r, i, j, k, b], where r = (i,j,k).transpose() * W * (i,j,k).

Definition at line 102 of file NeighborList.hh.

EIGEN_MAKE_ALIGNED_OPERATOR_NEW typedef Index CASM::PrimNeighborList::Scalar

Definition at line 30 of file NeighborList.hh.


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