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

#include <Structure.hh>

+ Inheritance diagram for CASM::Structure:

Detailed Description

Structure specifies the lattice and atomic basis of a crystal.

Definition at line 29 of file Structure.hh.

Public Member Functions

 Structure ()
 
 Structure (const Lattice &init_lat)
 
 Structure (const BasicStructure< Site > &base)
 
 Structure (const fs::path &filepath)
 
 Structure (const Structure &RHS)
 
const MasterSymGroupfactor_group () const
 
const SymGrouppoint_group () const
 
SymGroupRep const * basis_permutation_symrep () const
 
SymGroupRepID basis_permutation_symrep_ID () const
 
std::vector< Specieget_struc_specie () const
 Returns an Array of each possible Specie in this Structure. More...
 
std::vector< Moleculeget_struc_molecule () const
 Returns an Array of each possible Molecule in this Structure. More...
 
std::vector< std::string > get_struc_molecule_name () const
 Returns an Array of each possible Molecule in this Structure. More...
 
Eigen::VectorXi get_num_each_specie () const
 
Eigen::VectorXi get_num_each_molecule () const
 
Structureoperator= (const Structure &RHS)
 
void copy_attributes_from (const Structure &RHS)
 copy all non-derived members More...
 
void reset ()
 clears symmetry, site internals, and other attributes More...
 
void set_lattice (const Lattice &lattice, COORD_TYPE mode)
 
void generate_factor_group (double map_tol=TOL) const
 
void generate_factor_group_slow (double map_tol=TOL) const
 Obtain factor group by testing all operations of the lattice point_group and keep. More...
 
void fg_converge (double large_tol)
 generate factor groups for a range of tol values, prints results to screen (for now) More...
 
void fg_converge (double small_tol, double large_tol, double increment)
 
SymGroupRepID generate_basis_permutation_representation (bool verbose=false) const
 Obtain the basis permutation representation of factor_group, returns its rep_id, and sets internal basis_perm_rep_ID. More...
 
void symmetrize (const SymGroup &relaxed_factors)
 
void symmetrize (const double &tolerace)
 
void fill_supercell (const Structure &prim, double map_tol=TOL)
 fill an empty structure with the basis of its corresponding primitive cell - performs optimized factor_group expansion More...
 
Structure create_superstruc (const Lattice &scel_lat, double map_tol=TOL) const
 Shortcut routine to create a supercell structure and fill it with sites. More...
 
void map_superstruc_to_prim (Structure &prim)
 Figures out which prim basis each superstructure basis corresponds to. More...
 
void set_occs (Array< int > occ_index)
 Setting the current occupants of the structure to those specified by an array of integers. More...
 
void generate_basis_bouquet (const SiteOrbitree &in_tree, SiteOrbitree &out_tree, Index num_sites, double tol)
 
void generate_asym_bouquet (const SiteOrbitree &in_tree, SiteOrbitree &out_tree, Index num_sites, double tol)
 
void generate_flowertrees_safe (const SiteOrbitree &in_tree, Array< SiteOrbitree > &out_trees, double tol=TOL)
 Gets clusters of every size radiating from one site and saves them to a flowertree. A garland for each site is constructed. More...
 
void generate_flowertrees (const SiteOrbitree &in_tree, Array< SiteOrbitree > &out_trees, double tol=TOL)
 
Structure stack_on (const Structure &understruc, bool override=0) const
 
Structure get_reflection () const
 Return reflection of structure. More...
 
void clump_atoms (double maxdist, double tol)
 If atoms are too close together, average their distance and make them one. More...
 
void sort_basis ()
 Rearrange basis by grouping atoms by type. More...
 
Structure stamp_with (SiteCluster stamp, bool lat_override=0, bool im_override=0) const
 
Array< Structurebedazzle (Array< SiteCluster > stamps, bool lat_override=0, bool im_override=0) const
 
Array< Array< Array< double > > > get_NN_table (const double &maxr, SiteOrbitree &bouquet, double tol)
 
Array< Array< Array< double > > > get_NN_table (const double &maxr, double tol)
 
void add_vacuum_shift (Structure &new_surface_struc, double vacuum_thickness, Eigen::Vector3d shift, COORD_TYPE mode) const
 Add vacuum and shift c vector. The vacuum is always added parallel to c, and the shift vector should also be parallel to the ab plane (x,y,0) More...
 
void add_vacuum_shift (Structure &new_surface_struc, double vacuum_thickness, Coordinate shift) const
 
void add_vacuum (Structure &new_surface_struc, double vacuum_thickness) const
 Adds vacuum layer on top of ab plane. More...
 
Structureoperator+= (const Coordinate &shift)
 Translates all atoms in cell. More...
 
Structureoperator-= (const Coordinate &shift)
 
void intpol (Structure end_struc, int Nofimag, PERIODICITY_TYPE mode, Array< Structure > &images)
 
void print_site_symmetry (std::ostream &stream, COORD_TYPE mode, int shorttag, double tol)
 For each symmetrically distinct site, print the symmetry operations that map it onto itself. More...
 
bool read_species ()
 
void assign_species (Array< std::string > &names, Array< double > &masses, Array< double > &magmoms, Array< double > &Us, Array< double > &Js)
 
jsonParserto_json (jsonParser &json) const
 
void from_json (const jsonParser &json)
 
Index find (const CoordType2 &test_site, double tol=TOL) const
 
Index find (const CoordType2 &test_site, const Coordinate &shift, double tol) const
 
const Latticelattice () const
 
UnitCellCoord get_unit_cell_coord (const CoordType2 &test_site, double tol=TOL) const
 
void copy_attributes_from (const BasicStructure &RHS)
 
void update ()
 
void set_site_internals ()
 Associate each site with its basis index by setting its internal flags (asym_ind -> -1) More...
 
void within ()
 Translate all basis sites so that they are inside the unit cell. More...
 
Site get_site (const UnitCellCoord &ucc) const
 
void set_basis (Array< Site > basis_in)
 Manually set the basis sites. More...
 
void generate_factor_group (SymGroup &factor_group, double map_tol) const
 apply a symmetry operation to the current structure (may change lattice vectors and order of basis atoms) More...
 
void generate_factor_group_slow (SymGroup &factor_group, double map_tol) const
 
void fg_converge (SymGroup &factor_group, double small_tol, double large_tol, double increment)
 
SymGroupRepID generate_basis_permutation_representation (const MasterSymGroup &factor_group, bool verbose) const
 
bool is_primitive (double prim_tol=TOL) const
 
bool is_primitive (BasicStructure &new_prim, double prim_tol=TOL) const
 
void fill_supercell (const BasicStructure &prim, double map_tol=TOL)
 fill an empty structure with the basis of its corresponding primitive cell More...
 
void generate_flowertrees_safe (const SiteOrbitree &in_tree, Array< SiteOrbitree > &out_trees)
 Gets clusters of every size radiating from one site and saves them to a flowertree. A garland for each site is constructed. More...
 
void generate_flowertrees (const SiteOrbitree &in_tree, Array< SiteOrbitree > &out_trees)
 
void map_superstruc_to_prim (BasicStructure &prim, const SymGroup &point_group)
 Figures out which prim basis each superstructure basis corresponds to. More...
 
void merge_sites (double maxdist)
 If atoms are too close together, average their distance and make them one. More...
 
void add_vacuum_shift (BasicStructure &new_surface_struc, double vacuum_thickness, Eigen::Vector3d shift, COORD_TYPE mode) const
 Add vacuum and shift c vector. The vacuum is always added parallel to c, and the shift vector should also be parallel to the ab plane (x,y,0) More...
 
void add_vacuum_shift (BasicStructure &new_surface_struc, double vacuum_thickness, Coordinate shift) const
 
void add_vacuum (BasicStructure &new_surface_struc, double vacuum_thickness) const
 Adds vacuum layer on top of ab plane. More...
 
Index max_possible_vacancies () const
 Counts sites that allow vacancies. More...
 
const Latticeget_lattice () const
 Return a reference to the lattice. More...
 
virtual void read (std::istream &stream)
 Print intpolated images in seperate directries. More...
 
void print_xyz (std::ostream &stream) const
 Output other formats. More...
 
void print_cif (std::ostream &stream) const
 

Public Attributes

std::string title
 User-specified name of this Structure. More...
 
Array< Sitebasis
 Lattice vectors that specifies periodicity of the crystal. More...
 

Protected Attributes

MasterSymGroup m_factor_group
 
SymGroupRepID basis_perm_rep_ID
 This holds the representation id of the permutation representation. More...
 
bool SD_flag
 Specifies whether selectice dynamics is on or of for DFT calculations. More...
 
Lattice m_lattice
 

Private Member Functions

void main_print (std::ostream &stream, COORD_TYPE mode, bool version5, int option) const
 *** Inherited from BasicStructure<Site> More...
 

Constructor & Destructor Documentation

CASM::Structure::Structure ( )
inline

Definition at line 57 of file Structure.hh.

CASM::Structure::Structure ( const Lattice init_lat)
inlineexplicit

Definition at line 58 of file Structure.hh.

CASM::Structure::Structure ( const BasicStructure< Site > &  base)
inlineexplicit

Definition at line 59 of file Structure.hh.

CASM::Structure::Structure ( const fs::path &  filepath)
explicit

Definition at line 18 of file Structure.cc.

CASM::Structure::Structure ( const Structure RHS)

Have to explicitly define the copy constructor so that factor_group does not depend on the lattice of 'RHS'

Definition at line 31 of file Structure.cc.

Member Function Documentation

void CASM::BasicStructure< Site >::add_vacuum ( BasicStructure< Site > &  new_surface_struc,
double  vacuum_thickness 
) const
inherited

Adds vacuum layer on top of ab plane.

void CASM::Structure::add_vacuum ( Structure new_surface_struc,
double  vacuum_thickness 
) const

Adds vacuum layer on top of ab plane.

Definition at line 1117 of file Structure.cc.

void CASM::BasicStructure< Site >::add_vacuum_shift ( BasicStructure< Site > &  new_surface_struc,
double  vacuum_thickness,
Eigen::Vector3d  shift,
COORD_TYPE  mode 
) const
inherited

Add vacuum and shift c vector. The vacuum is always added parallel to c, and the shift vector should also be parallel to the ab plane (x,y,0)

Call this on a structure to get new_surface_struc: the structure with a layer of vacuum added parallel to the ab plane. vacuum_thickness: thickness of vacuum layer (Angstroms) shift: shift vector from layer to layer, assumes FRAC unless specified. The shift vector should only have values relative to a and b vectors (eg x, y, 0). Default shift is zero.

void CASM::BasicStructure< Site >::add_vacuum_shift ( BasicStructure< Site > &  new_surface_struc,
double  vacuum_thickness,
Coordinate  shift 
) const
inherited
void CASM::Structure::add_vacuum_shift ( Structure new_surface_struc,
double  vacuum_thickness,
Eigen::Vector3d  shift,
COORD_TYPE  mode 
) const

Add vacuum and shift c vector. The vacuum is always added parallel to c, and the shift vector should also be parallel to the ab plane (x,y,0)

Call this on a structure to get new_surface_struc: the structure with a layer of vacuum added parallel to the ab plane. vacuum_thickness: thickness of vacuum layer (Angstroms) shift: shift vector from layer to layer, assumes FRAC unless specified. The shift vector should only have values relative to a and b vectors (eg x, y, 0). Default shift is zero.

Definition at line 1088 of file Structure.cc.

void CASM::Structure::add_vacuum_shift ( Structure new_surface_struc,
double  vacuum_thickness,
Coordinate  shift 
) const

Definition at line 1110 of file Structure.cc.

void CASM::Structure::assign_species ( Array< std::string > &  names,
Array< double > &  masses,
Array< double > &  magmoms,
Array< double > &  Us,
Array< double > &  Js 
)
SymGroupRep const * CASM::Structure::basis_permutation_symrep ( ) const

Definition at line 107 of file Structure.cc.

SymGroupRepID CASM::Structure::basis_permutation_symrep_ID ( ) const

Definition at line 113 of file Structure.cc.

Array< Structure > CASM::Structure::bedazzle ( Array< SiteCluster stamps,
bool  lat_override = 0,
bool  im_override = 0 
) const

Same as stamp_with, but takes an array of clusters and returns and array of structures. Basically saves you a for loop. Perhaps this should actually take an OrbitBranch?

Definition at line 945 of file Structure.cc.

void CASM::Structure::clump_atoms ( double  mindist,
double  tol 
)

If atoms are too close together, average their distance and make them one.

Given a minimum distance, find pair clusters within this length and average out the distance between the two atoms and put one atom there, eliminating the other two. Meant to be used for grain boundaries.

This function is meant for averaging atoms of the same type!

Definition at line 793 of file Structure.cc.

void CASM::Structure::copy_attributes_from ( const Structure RHS)

copy all non-derived members

Definition at line 62 of file Structure.cc.

void CASM::BasicStructure< Site >::copy_attributes_from ( const BasicStructure< Site > &  RHS)
inherited

Use this is the copy interface for things that derive from BasicStructure<> It should be overloaded in derived classes so that all important attributes besides lattice, basis, and title get copied

Structure CASM::Structure::create_superstruc ( const Lattice scel_lat,
double  map_tol = TOL 
) const

Shortcut routine to create a supercell structure and fill it with sites.

Operates on the primitive structure and takes as an argument the supercell lattice. It then returns a new superstructure.

This is similar to the Lattice*Primitive routine which returns a new superstructure. Unlike the fill_supercell routine which takes the primitive structure, this WILL fill the sites.

Definition at line 335 of file Structure.cc.

const MasterSymGroup & CASM::Structure::factor_group ( ) const

Definition at line 94 of file Structure.cc.

void CASM::Structure::fg_converge ( double  large_tol)

generate factor groups for a range of tol values, prints results to screen (for now)

Definition at line 231 of file Structure.cc.

void CASM::Structure::fg_converge ( double  small_tol,
double  large_tol,
double  increment 
)

Definition at line 225 of file Structure.cc.

void CASM::BasicStructure< Site >::fg_converge ( SymGroup factor_group,
double  small_tol,
double  large_tol,
double  increment 
)
inherited
void CASM::Structure::fill_supercell ( const Structure prim,
double  map_tol = TOL 
)

fill an empty structure with the basis of its corresponding primitive cell - performs optimized factor_group expansion

Returns true if the structure describes a crystal primitive cell i.e., no translation smaller than a lattice vector can map the structure onto itself Returns true if the structure describes a crystal primitive cell and finds the primitive cell and stores it in 'new_prim'

It is NOT wise to use this function unless you have already initialized a superstructure with lattice vectors.

It is more wise to use the two methods that call this method: Either the overloaded * operator which does: SCEL_Lattice * Prim_Structrue = New_Superstructure — or — New_Superstructure=Prim_Structure.create_superstruc(SCEL_Lattice);

Both of these will return NEW superstructures.

Definition at line 263 of file Structure.cc.

void CASM::BasicStructure< Site >::fill_supercell ( const BasicStructure< Site > &  prim,
double  map_tol = TOL 
)
inherited

fill an empty structure with the basis of its corresponding primitive cell

It is NOT wise to use this function unless you have already initialized a superstructure with lattice vectors.

It is more wise to use the two methods that call this method: Either the overloaded * operator which does: SCEL_Lattice * Prim_Structrue = New_Superstructure — or — New_Superstructure=Prim_BasicStructure<CoordType>.create_superstruc(SCEL_Lattice);

Both of these will return NEW superstructures.

Index CASM::BasicStructure< Site >::find ( const CoordType2 &  test_site,
double  tol = TOL 
) const
inherited

return basis index of site that matches test_coord, if it is in basis otherwise, returns basis.size()

Index CASM::BasicStructure< Site >::find ( const CoordType2 &  test_site,
const Coordinate shift,
double  tol 
) const
inherited

return basis index of site that matches test_site+shift, if it is in basis otherwise, returns basis.size()

void CASM::Structure::from_json ( const jsonParser json)

re-throw exceptions

Definition at line 1241 of file Structure.cc.

void CASM::Structure::generate_asym_bouquet ( const SiteOrbitree in_tree,
SiteOrbitree out_tree,
Index  num_sites,
double  tol 
)

For each asymmetric unit site, find each cluster in 'in_tree' that contain that basis site, create an orbit of all equivalent clusters that also contain the basis site and store each of these orbits in out_tree num_sites specifies what size of clusters are desired (e.g., pairs, triplets, etc)

Fills out_tree with OrbitBranches with pivots corresponding to the different asym_unit sites. The orbits of each branch are "petals" radiating from that branch's pivot.

Parameters
in_treeGeneral tree that holds clusters of all sizes
out_treeHolds all flowers (orbits) radiating from asym_unit sites
num_sitesSize of clusters filling up out_tree

Definition at line 501 of file Structure.cc.

void CASM::Structure::generate_basis_bouquet ( const SiteOrbitree in_tree,
SiteOrbitree out_tree,
Index  num_sites,
double  tol 
)

For each basis site, find each cluster in 'in_tree' that contain that basis site, create an orbit of all equivalent clusters that also contain the basis site and store each of these orbits in out_tree num_sites specifies what size of clusters are desired (e.g., pairs, triplets, etc)

Definition at line 446 of file Structure.cc.

SymGroupRepID CASM::Structure::generate_basis_permutation_representation ( bool  verbose = false) const

Obtain the basis permutation representation of factor_group, returns its rep_id, and sets internal basis_perm_rep_ID.

Definition at line 1182 of file Structure.cc.

SymGroupRepID CASM::BasicStructure< Site >::generate_basis_permutation_representation ( const MasterSymGroup factor_group,
bool  verbose 
) const
inherited
void CASM::Structure::generate_factor_group ( double  map_tol = TOL) const

apply a symmetry operation to the current structure (does not change lattice vectors, because 'op' is assumed to be a lattice homomorphism) determines primitive cell, finds its factor group using generate_factor_group_slow, and then expands the factor group into the supercell using lattice translations

Definition at line 85 of file Structure.cc.

void CASM::BasicStructure< Site >::generate_factor_group ( SymGroup factor_group,
double  map_tol 
) const
inherited

apply a symmetry operation to the current structure (may change lattice vectors and order of basis atoms)

void CASM::Structure::generate_factor_group_slow ( double  map_tol = TOL) const

Obtain factor group by testing all operations of the lattice point_group and keep.

Definition at line 77 of file Structure.cc.

void CASM::BasicStructure< Site >::generate_factor_group_slow ( SymGroup factor_group,
double  map_tol 
) const
inherited
void CASM::BasicStructure< Site >::generate_flowertrees ( const SiteOrbitree in_tree,
Array< SiteOrbitree > &  out_trees 
)
inherited
void CASM::Structure::generate_flowertrees ( const SiteOrbitree in_tree,
Array< SiteOrbitree > &  out_trees,
double  tol = TOL 
)

A flowertree has clusters of every size organized like the branches of a regular orbitree, but every cluster in the entire flowertree has the same pivot.

Definition at line 410 of file Structure.cc.

void CASM::BasicStructure< Site >::generate_flowertrees_safe ( const SiteOrbitree in_tree,
Array< SiteOrbitree > &  out_trees 
)
inherited

Gets clusters of every size radiating from one site and saves them to a flowertree. A garland for each site is constructed.

void CASM::Structure::generate_flowertrees_safe ( const SiteOrbitree in_tree,
Array< SiteOrbitree > &  out_trees,
double  tol = TOL 
)

Gets clusters of every size radiating from one site and saves them to a flowertree. A garland for each site is constructed.

A flowertree is NOT a bouquet. Bouquets are constructed from an orbitree (generate_basis_bouquet) by taking every cluster of a given size and storing them in branches, where each branch contains clusters with a common pivot. A flowertree has clusters of every size organized like the branches of a regular orbitree, but every cluster in the entire flowertree has the same pivot.

Begin by making n bouquets, where n is the size of the larges cluster in the provided orbitree. Then shuffle all the bouquets together and sort them by basis site. The result is an array of "flowertrees", which have been made by picking out flowers of every bouquet that have the same pivot. You end up with one flowertree per basis site each of which has clusters from size 1 to n.

Definition at line 562 of file Structure.cc.

const Lattice& CASM::BasicStructure< Site >::get_lattice ( ) const
inlineinherited

Return a reference to the lattice.

Definition at line 177 of file BasicStructure.hh.

Array< Array< Array< double > > > CASM::Structure::get_NN_table ( const double &  maxr,
SiteOrbitree bouquet,
double  tol 
)

Goes to a specified site of the basis and makes a flower tree of pairs. It then stores the length and multiplicity of the pairs in a double array, giving you a strict nearest neighbor table. This version also fills up a SiteOrbitree in case you want to keep it. Blatantly copied from Anna's routine in old new CASM

Definition at line 965 of file Structure.cc.

Array< Array< Array< double > > > CASM::Structure::get_NN_table ( const double &  maxr,
double  tol 
)

Goes to a specified site of the basis and makes a flower tree of pairs. It then stores the length and multiplicity of the pairs in a double array, giving you a strict nearest neighbor table. The bouquet used for this falls into the void.

Definition at line 1008 of file Structure.cc.

Eigen::VectorXi CASM::Structure::get_num_each_molecule ( ) const

Returns a list of how many of each molecule exist in this Structure The molecule types are ordered according to get_struc_molecule()

Definition at line 206 of file Structure.cc.

Eigen::VectorXi CASM::Structure::get_num_each_specie ( ) const

Returns a list of how many of each specie exist in this Structure The Specie types are ordered according to get_struc_specie()

Definition at line 184 of file Structure.cc.

Structure CASM::Structure::get_reflection ( ) const

Return reflection of structure.

Makes reflection of structure and returns new reflected strucutre. Meant for primitive structures, but can be forced on superstructures. If forced on non primitive structure then the reflected structure will think it's primitive. To avoid this issue reflect primitive structure and then make a reflected superlattice to fill it with.

Definition at line 769 of file Structure.cc.

Site CASM::BasicStructure< Site >::get_site ( const UnitCellCoord ucc) const
inherited
std::vector< Molecule > CASM::Structure::get_struc_molecule ( ) const

Returns an Array of each possible Molecule in this Structure.

Definition at line 146 of file Structure.cc.

std::vector< std::string > CASM::Structure::get_struc_molecule_name ( ) const

Returns an Array of each possible Molecule in this Structure.

Definition at line 166 of file Structure.cc.

std::vector< Specie > CASM::Structure::get_struc_specie ( ) const

Returns an Array of each possible Specie in this Structure.

Definition at line 123 of file Structure.cc.

UnitCellCoord CASM::BasicStructure< Site >::get_unit_cell_coord ( const CoordType2 &  test_site,
double  tol = TOL 
) const
inherited

Return the UnitCellCoord corresponding to test_site (i.e., finds the basis index and the lattice translation)

Using the lattice of (*this), this function will return a UnitCellCoord that corresponds to a site passed to it within a given tolerance. This function is useful for making a nearest neighbor table from sites that land outside of the primitive cell.

void CASM::Structure::intpol ( Structure  end_struc,
int  Nofimag,
PERIODICITY_TYPE  mode,
Array< Structure > &  images 
)

Creates Nofimag POSCAR files by interpolating linearly between structures in star_stru and end_stru for exact interpolation, choose "LOCAL" or "1", for nearest-image interpolation, choose "PERIODIC" or "0"

Creates N of image POSCAR files in seperate directories by interpolating linearly between structures. For exact interpolation, choose " LOCAL " or "1" in mode For nearest-image interpolation, chosse "PERIODIC" or "0"

Stored in Array<Structure> images; images[0] = start structure, images[Nofimage+1] = end structure

Definition at line 1139 of file Structure.cc.

bool CASM::BasicStructure< Site >::is_primitive ( double  prim_tol = TOL) const
inherited

Returns true if the structure describes a crystal primitive cell i.e., no translation smaller than a lattice vector can map the structure onto itself

Determines if structure is primitive description of the crystal

bool CASM::BasicStructure< Site >::is_primitive ( BasicStructure< Site > &  new_prim,
double  prim_tol = TOL 
) const
inherited

Returns true if the structure describes a crystal primitive cell and finds the primitive cell and stores it in 'new_prim'

Determines if structure is primitive description of the crystal If not, finds primitive cell and copies to new_prim

const Lattice& CASM::BasicStructure< Site >::lattice ( ) const
inlineinherited

Definition at line 74 of file BasicStructure.hh.

void CASM::Structure::main_print ( std::ostream &  stream,
COORD_TYPE  mode,
bool  version5,
int  option 
) const
private

*** Inherited from BasicStructure<Site>

void CASM::Structure::map_superstruc_to_prim ( Structure prim)

Figures out which prim basis each superstructure basis corresponds to.

Definition at line 605 of file Structure.cc.

void CASM::BasicStructure< Site >::map_superstruc_to_prim ( BasicStructure< Site > &  prim,
const SymGroup point_group 
)
inherited

Figures out which prim basis each superstructure basis corresponds to.

Index CASM::BasicStructure< Site >::max_possible_vacancies ( ) const
inherited

Counts sites that allow vacancies.

void CASM::BasicStructure< Site >::merge_sites ( double  maxdist)
inherited

If atoms are too close together, average their distance and make them one.

Structure & CASM::Structure::operator+= ( const Coordinate shift)

Translates all atoms in cell.

Definition at line 1201 of file Structure.cc.

Structure & CASM::Structure::operator-= ( const Coordinate shift)

Definition at line 1213 of file Structure.cc.

Structure & CASM::Structure::operator= ( const Structure RHS)

Have to explicitly define the assignment operator so that sites in this structure do not depend on the lattice of 'RHS'

Definition at line 40 of file Structure.cc.

const SymGroup & CASM::Structure::point_group ( ) const

Definition at line 101 of file Structure.cc.

void CASM::BasicStructure< Site >::print_cif ( std::ostream &  stream) const
inherited
void CASM::Structure::print_site_symmetry ( std::ostream &  stream,
COORD_TYPE  mode,
int  shorttag = 1,
double  tol = TOL 
)

For each symmetrically distinct site, print the symmetry operations that map it onto itself.

Definition at line 348 of file Structure.cc.

void CASM::BasicStructure< Site >::print_xyz ( std::ostream &  stream) const
inherited

Output other formats.

virtual void CASM::BasicStructure< Site >::read ( std::istream &  stream)
virtualinherited

Print intpolated images in seperate directries.

bool CASM::Structure::read_species ( )
void CASM::Structure::reset ( )
virtual

clears symmetry, site internals, and other attributes

Should we also invalidate the occupants? for(Index i = 0; i < basis.size(); i++) { basis[i].site_occupant.set_value(-1); }

Reimplemented from CASM::BasicStructure< Site >.

Definition at line 385 of file Structure.cc.

void CASM::BasicStructure< Site >::set_basis ( Array< Site basis_in)
inherited

Manually set the basis sites.

Allows for the basis elements of a basic structure to be manually set, e.g. as in jsonParser.cc.

void CASM::Structure::set_lattice ( const Lattice lattice,
COORD_TYPE  mode 
)

change the lattice and update site coordinates. Argument 'mode' specifies which mode is preserved e.g.: struc.set_lattice(new_lat, CART) calculates all Cartesian coordinates, invalidates the FRAC coordinates, and changes the lattice

Definition at line 664 of file Structure.cc.

void CASM::Structure::set_occs ( Array< int >  occ_index)

Setting the current occupants of the structure to those specified by an array of integers.

Definition at line 1190 of file Structure.cc.

void CASM::BasicStructure< Site >::set_site_internals ( )
inherited

Associate each site with its basis index by setting its internal flags (asym_ind -> -1)

void CASM::Structure::sort_basis ( )

Rearrange basis by grouping atoms by type.

Loop through basis and rearrange atoms by type. Uses bubble sort algorithm by comparing integer values of the strings assigned to the basis occupants.

The basis gets sorted in a sort of alphabetical way, so be mindful of any POTCARs you might have if you run this.

Definition at line 858 of file Structure.cc.

Structure CASM::Structure::stack_on ( const Structure understruc,
bool  override = 0 
) const

Returns a heterostructre. The structure this function is called onto will be stretched to match a and b of the structure is it passed (understruc). The returned strucure will contain the basis of both structures, stacked on top of each other in the c direction.

This function only matches the ab planes.

Definition at line 692 of file Structure.cc.

Structure CASM::Structure::stamp_with ( SiteCluster  stamp,
bool  lat_override = 0,
bool  im_override = 0 
) const

Decorate your structure easily! Given a cluster, this routine will copy your structure and replace sites of the copy with sites from the cluster. Bedazzle! The lattice of your structure and cluster must be related by an integer matrix (supercell). Be mindful with the size of your stamped structure, it has to be big enough to fit the cluster.

Definition at line 888 of file Structure.cc.

void CASM::Structure::symmetrize ( const SymGroup relaxed_factors)

Given a symmetry group, the basis of the structure will have each operation applied to it. The resulting set of basis from performing these operations will be averaged out, yielding a new average basis that will replace the current one.

Definition at line 1022 of file Structure.cc.

void CASM::Structure::symmetrize ( const double &  tolerance)

Same as the other symmetrize routine, except this one assumes that the symmetry group you mean to use is the factor group of your structure within a certain tolerance.

Notice that the tolerance is also used on your point group!!

Definition at line 1069 of file Structure.cc.

jsonParser & CASM::Structure::to_json ( jsonParser json) const

Definition at line 1224 of file Structure.cc.

void CASM::BasicStructure< Site >::update ( )
inherited
void CASM::BasicStructure< Site >::within ( )
inherited

Translate all basis sites so that they are inside the unit cell.

Member Data Documentation

Array<Site > CASM::BasicStructure< Site >::basis
inherited

Lattice vectors that specifies periodicity of the crystal.

Definition at line 42 of file BasicStructure.hh.

SymGroupRepID CASM::Structure::basis_perm_rep_ID
mutableprotected

This holds the representation id of the permutation representation.

Definition at line 36 of file Structure.hh.

MasterSymGroup CASM::Structure::m_factor_group
mutableprotected

Group symmetry operations that map the lattice and basis of Structure onto themselves, assuming that the crystal is periodic

Definition at line 34 of file Structure.hh.

Lattice CASM::BasicStructure< Site >::m_lattice
protectedinherited

Definition at line 32 of file BasicStructure.hh.

bool CASM::Structure::SD_flag
protected

Specifies whether selectice dynamics is on or of for DFT calculations.

Definition at line 38 of file Structure.hh.

std::string CASM::BasicStructure< Site >::title
inherited

User-specified name of this Structure.

Definition at line 39 of file BasicStructure.hh.


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