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

#include <NeighborList.hh>

Detailed Description

The SuperNeighborList gives the linear indices of neighboring sites and unitcells in a particular Supercell.

  • The SuperNeighborList is constructed from a PrimGrid and PrimNeighborList
  • The linear index of the neighboring unit cell is determined for all unit cells in a supercell
  • The linear index of sites in each neighboring unit cell is determined for all unit cells in a supercell

Definition at line 125 of file NeighborList.hh.

Public Types

typedef Index size_type
 

Public Member Functions

 SuperNeighborList (const PrimGrid &prim_grid, const PrimNeighborList &prim_nlist)
 Constructor. More...
 
size_type unitcell_index (size_type site_index) const
 get unit cell index from site index More...
 
size_type sublat_index (size_type site_index) const
 get sublat 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
< SuperNeighborList
clone () 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

Definition at line 129 of file NeighborList.hh.

Constructor & Destructor Documentation

CASM::SuperNeighborList::SuperNeighborList ( const PrimGrid prim_grid,
const PrimNeighborList prim_nlist 
)

Constructor.

Parameters
prim_gridA grid of unit cells describing the supercell this neighbor list pertains to
prim_nlist_begin,prim_nlist_endIterators over range [begin, end) of UnitCell that are the neighbors of the origin UnitCell
sublat_begin,sublat_endIterators over range [begin, end) of basis sites that should be included as UnitCellCoord neighbors
  • 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 189 of file NeighborList.cc.

Member Function Documentation

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

Clone.

Definition at line 247 of file NeighborList.cc.

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 236 of file NeighborList.cc.

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 224 of file NeighborList.cc.

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

get sublat index from site index

Definition at line 140 of file NeighborList.hh.

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

get unit cell index from site index

Definition at line 135 of file NeighborList.hh.

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 229 of file NeighborList.cc.

Member Data Documentation

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 160 of file NeighborList.hh.

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 + n*scel.volume() are identical.
  • So m_site.size() == m_prim_grid_size

Definition at line 168 of file NeighborList.hh.

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

m_unitcell[unitcell_index][neighbor unitcell index]

Definition at line 171 of file NeighborList.hh.


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