CASM  1.1.0
A Clusters Approach to Statistical Mechanics
CASM::SuperNeighborList Class Reference

#include <NeighborList.hh>

Detailed Description

SuperNeighborList, linear indices of neighboring sites and unit cells

The SuperNeighborList takes the ordering of unit cells neighboring the origin unit cell from the PrimNeighborList and uses it to construct lists of linear indices of neighboring unit cells and neighboring sites for each unit cell in a particular supercell.

Neighbor list index naming conventions:

Definition at line 177 of file NeighborList.hh.

Public Types

typedef Index size_type
 

Public Member Functions

 SuperNeighborList (Eigen::Matrix3l const &transformation_matrix_to_super, PrimNeighborList const &prim_nlist)
 Constructor. More...
 
 SuperNeighborList (const xtal::Superlattice &prim_grid, const PrimNeighborList &prim_nlist)
 Constructor. More...
 
size_type unitcell_index (size_type site_index) const
 Get unitcell_index from site_index. More...
 
size_type sublat_index (size_type site_index) const
 Get sublattice index from site_index. More...
 
const std::vector< size_type > & sites (size_type unitcell_index) const
 const Access the list of sites neighboring a particular unit cell More...
 
const std::vector< size_type > & unitcells (size_type unitcell_index) const
 const Access the list of unitcells neighboring a particular unit cell More...
 
bool overlaps () const
 Returns true if periodic images of the neighbor list overlap. More...
 
std::unique_ptr< SuperNeighborListclone () const
 Clone. More...
 

Private Attributes

size_type m_prim_grid_size
 store prim grid size for site index -> unitcell index conversion unitcell_index = site_index % m_prim_grid_size More...
 
std::vector< std::vector< size_type > > m_site
 m_site[unitcell_index][neighbor site index] More...
 
std::vector< std::vector< size_type > > m_unitcell
 m_unitcell[unitcell_index][neighbor unitcell index] More...
 

Member Typedef Documentation

◆ size_type

Definition at line 179 of file NeighborList.hh.

Constructor & Destructor Documentation

◆ SuperNeighborList() [1/2]

CASM::SuperNeighborList::SuperNeighborList ( Eigen::Matrix3l const &  transformation_matrix_to_super,
PrimNeighborList const &  prim_nlist 
)

Constructor.

Parameters
transformation_matrix_to_superThe transformation matrix, T, from primitive lattice vectors to supercell lattice vectors: supercell_lat_column_mat = T * prim_lat_column_mat.
prim_nlistA reference to a PrimNeighborList defining the neighborhood for the origin unit cell. It is only used in the constructor, any changes to the PrimNeighborList will not be reflected in this SuperNeighborList.
  • 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]
  • The sublattice iterators enable restricting the neighbor list to only sites that have degrees of freedom

Definition at line 235 of file NeighborList.cc.

◆ SuperNeighborList() [2/2]

CASM::SuperNeighborList::SuperNeighborList ( const xtal::Superlattice superlattice,
const PrimNeighborList prim_nlist 
)

Constructor.

Parameters
superlatticeThe superlattice this neighbor list pertains to, coupled with its primitive tiling unit
prim_nlistA reference to a PrimNeighborList defining the neighborhood for the origin unit cell. It is only used in the constructor, any changes to the PrimNeighborList will not be reflected in this SuperNeighborList.
  • 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]
  • The sublattice iterators enable restricting the neighbor list to only sites that have degrees of freedom

Definition at line 289 of file NeighborList.cc.

Member Function Documentation

◆ clone()

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

Clone.

Definition at line 320 of file NeighborList.cc.

◆ overlaps()

bool CASM::SuperNeighborList::overlaps ( ) const

Returns true if periodic images of the neighbor list overlap.

If periodic images of the neighborhood overlap, Clexulator 'delta' values will be incorrect.

Definition at line 310 of file NeighborList.cc.

◆ sites()

const std::vector< SuperNeighborList::size_type > & CASM::SuperNeighborList::sites ( size_type  unitcell_index) const

const Access the list of sites neighboring a particular unit cell

const Access the list of sites neighboring a particular unitcell

Definition at line 295 of file NeighborList.cc.

◆ sublat_index()

size_type CASM::SuperNeighborList::sublat_index ( size_type  site_index) const
inline

Get sublattice index from site_index.

Definition at line 195 of file NeighborList.hh.

◆ unitcell_index()

size_type CASM::SuperNeighborList::unitcell_index ( size_type  site_index) const
inline

Get unitcell_index from site_index.

Definition at line 190 of file NeighborList.hh.

◆ unitcells()

const std::vector< SuperNeighborList::size_type > & CASM::SuperNeighborList::unitcells ( size_type  unitcell_index) const

const Access the list of unitcells neighboring a particular unit cell

const Access the list of unitcells neighboring a particular unitcell

Definition at line 301 of file NeighborList.cc.

Member Data Documentation

◆ m_prim_grid_size

size_type CASM::SuperNeighborList::m_prim_grid_size
private

store prim grid size for site index -> unitcell index conversion unitcell_index = site_index % m_prim_grid_size

Definition at line 215 of file NeighborList.hh.

◆ m_site

std::vector<std::vector<size_type> > CASM::SuperNeighborList::m_site
private

m_site[unitcell_index][neighbor site index]

  • Configuration sites are ordered in blocks corresponding to each sublattice, b. Neighbors are specific only to a unitcell, so the neighbor list for site index s and s
  • So m_site.size() == m_prim_grid_size

Definition at line 225 of file NeighborList.hh.

◆ m_unitcell

std::vector<std::vector<size_type> > CASM::SuperNeighborList::m_unitcell
private

m_unitcell[unitcell_index][neighbor unitcell index]

Definition at line 228 of file NeighborList.hh.


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