CASM  1.1.0
A Clusters Approach to Statistical Mechanics
CASM::xtal Namespace Reference

Namespaces

 canonical
 
 Coordinate_impl
 
 HermiteCounter_impl
 
 impl
 
 Local
 
 SimpleStructureTools
 
 SpeciesAttribute_impl
 
 StrucMapping
 

Classes

class  BasicStructure
 BasicStructure specifies the lattice and atomic basis of a crystal. More...
 
class  Coordinate
 Represents cartesian and fractional coordinates. More...
 
class  COORD_MODE
 COORD_MODE specifies the current coordinate mode (Fractional or Cartesian) More...
 
class  PERIODICITY_MODE
 
class  DoFSet
 
class  SiteDoFSet
 
struct  DoFSetIsEquivalent_f
 
struct  SiteDoFSetIsEquivalent_f
 
class  HermiteCounter
 
class  IntegralCoordinateWithin_f
 
class  Lattice
 
class  LatticeIsEquivalent
 Lattice comparisons. More...
 
class  IsPointGroupOp
 Checks if operations are point group operations. More...
 
class  StrainCostCalculator
 
class  LatticeMap
 
class  UnitCellCoordIndexConverter
 
class  UnitCellIndexConverter
 
struct  _LinearIndexConverter
 
struct  _LinearIndexConverter< UnitCell >
 
struct  _LinearIndexConverter< UnitCellCoord >
 
class  AtomPosition
 An atomic species associated with a position in space. More...
 
class  Molecule
 Class representing a Molecule. More...
 
class  NiggliRep
 
class  StandardOrientationCompare
 
class  SimpleStrucMapCalculator
 
class  SimpleStructure
 Representation of a crystal of molecular and/or atomic occupants, and any additional properties. It does not require a primitive Structure or BasicStructure object to act as a reference and can be thought of as the class that embodies all the information of a config.json or properties.calc.json file. More...
 
class  Site
 
class  SpeciesAttribute
 
class  StrucMapCalculatorInterface
 
struct  LatticeNode
 Class describing the lattice-mapping portion of a particular mapping A general map for child_struc onto parent_struc may require forming a supercell of parent_struc (most commonly) and/or of child_struc. As such, the LatticeNode is specified in terms of superlattices of both parent_struc and child_struc, as well as deformation and rotation information sufficient to fully define the lattice map. More...
 
struct  AssignmentNode
 Structure to encode the solution of a constrained atomic assignmnet problem This describes the permutation, translation, and time-reversal of the atoms of a child structure to bring them into registration with the atoms of a parent structure (assuming periodic boundary conditions). Also records the constrained and unconstrained assignment costs. More...
 
struct  MappingNode
 
class  StrucMapper
 
class  Superlattice
 
class  ScelEnumProps
 Data structure for holding supercell enumeration properties. More...
 
class  SuperlatticeIterator
 Iterators used with SuperlatticeEnumerator. More...
 
class  SuperlatticeEnumerator
 A fake container of supercell matrices. More...
 
struct  SymOp
 
struct  SymOpCompare_f
 
struct  SymOpPeriodicCompare_f
 
struct  SymOpMatrixCompare_f
 
class  UnitCell
 Unit Cell Indices. More...
 
struct  Translatable
 CRTP class to implement '-=', '+', and '-' in terms of '+='. More...
 
class  UnitCellCoord
 Unit Cell Coordinates. More...
 

Typedefs

typedef std::vector< SymOpSymOpVector
 
template<typename IntegralCoordinateType >
using LinearIndexConverter = typename _LinearIndexConverter< IntegralCoordinateType >::type
 
typedef Eigen::Matrix3d SymOpMatrixType
 
typedef Eigen::Vector3d SymOpTranslationType
 
typedef bool SymOpTimeReversalType
 

Functions

IntegralCoordinateWithin_f make_bring_within_f (const Supercell &scel)
 Make IntegralCoordinateWithin_f for Supercell [deprecated]. More...
 
std::vector< UnitCellCoordsymop_site_map (SymOp const &_op, BasicStructure const &_struc)
 To which site a SymOp transforms each basis site. More...
 
template<typename ExternSymOp >
std::vector< UnitCellCoordsymop_site_map (ExternSymOp const &_op, BasicStructure const &_struc)
 To which site a SymOp transforms each basis site. More...
 
std::vector< UnitCellCoordsymop_site_map (SymOp const &_op, BasicStructure const &_struc, double _tol)
 To which site a SymOp transforms each basis site. More...
 
template<typename ExternSymOp >
std::vector< UnitCellCoordsymop_site_map (ExternSymOp const &_op, BasicStructure const &_struc, double _tol)
 To which site a SymOp transforms each basis site. More...
 
std::vector< Moleculestruc_molecule (BasicStructure const &_struc)
 Returns an Array of each possible Molecule in this Structure. More...
 
std::vector< std::string > struc_species (BasicStructure const &_struc)
 Returns an Array of each possible AtomSpecie in this Structure. More...
 
std::vector< std::string > struc_molecule_name (BasicStructure const &_struc)
 Returns an Array of each possible Molecule in this Structure. More...
 
std::vector< std::vector< std::string > > allowed_molecule_unique_names (BasicStructure const &_struc)
 Returns an Array of each possible Molecule in this Structure. More...
 
std::vector< std::vector< std::string > > allowed_molecule_names (BasicStructure const &_struc)
 Returns a vector with a list of allowed molecule names at each site. More...
 
std::vector< DoFKeyall_local_dof_types (BasicStructure const &_struc)
 
std::vector< DoFKeycontinuous_local_dof_types (BasicStructure const &_struc)
 
std::vector< DoFKeyglobal_dof_types (BasicStructure const &_struc)
 
std::vector< DoFKeyall_dof_types (BasicStructure const &_struc)
 
std::map< DoFKey, Indexlocal_dof_dims (BasicStructure const &_struc)
 
std::map< DoFKey, Indexglobal_dof_dims (BasicStructure const &_struc)
 
Index local_dof_dim (DoFKey const &_name, BasicStructure const &_struc)
 
bool has_strain_dof (BasicStructure const &structure)
 
DoFKey get_strain_dof_key (BasicStructure const &structure)
 
std::string get_strain_metric (DoFKey strain_dof_key)
 
Index find_index (const std::vector< Site > &basis, const Site &test_site, double tol)
 
bool is_primitive (const BasicStructure &struc, double tol=TOL)
 
BasicStructure make_primitive (const BasicStructure &non_primitive_struc, double tol=TOL)
 Returns the smallest possible tiling unit of the given structure. More...
 
std::pair< double, Eigen::Vector3d > calc_rotation_angle_and_axis (const SymOp &op, const Lattice &lat)
 Calculates the rotation angle and axis of a symmetry operation. This function is almost exactly identical to the constructor of SymInfo::SymInfo. More...
 
void sort_factor_group (std::vector< SymOp > &factor_group, const Lattice &lat)
 Sort SymOp in the SymGroup. More...
 
std::vector< SymOpmake_factor_group (const BasicStructure &struc, double tol=TOL)
 
std::vector< Eigen::PermutationMatrix< Eigen::Dynamic, Eigen::Dynamic, Index > > make_permutation_representation (const xtal::BasicStructure &struc, const std::vector< SymOp > &factor_group)
 Create the permutation group of a structure. More...
 
std::set< std::set< Index > > make_asymmetric_unit (const xtal::BasicStructure &struc, const std::vector< SymOp > &factor_group)
 Return indices of equivalent basis sites. More...
 
std::set< std::set< Index > > make_asymmetric_unit (const xtal::BasicStructure &struc)
 Return indices of equivalent basis sites. More...
 
BasicStructure symmetrize (const BasicStructure &structure, const std::vector< SymOp > &enforced_group)
 
template<typename IntegralType , int Options = 0>
BasicStructure make_superstructure (const BasicStructure &tiling_unit, const Eigen::Matrix< IntegralType, 3, 3, Options > &transformation_matrix)
 
Coordinate operator* (const SymOp &LHS, const Coordinate &RHS)
 
Coordinate operator+ (const Coordinate &LHS, const Coordinate &RHS)
 
Coordinate operator- (const Coordinate &LHS, const Coordinate &RHS)
 
std::vector< std::string > component_descriptions (DoFSet const &dofset)
 
template<typename DoFSetType >
std::map< std::string, DoFSetType > make_dofset_map (std::vector< DoFSetType > const &dofset_vec)
 
Eigen::MatrixXd dofset_transformation_matrix (const Eigen::MatrixXd &from_basis, const Eigen::MatrixXd &to_basis, double tol)
 Create the symmtery representation for going from one basis to another. More...
 
std::vector< UnitCellmake_lattice_points (const Eigen::Matrix3l &transformation_matrix)
 
std::vector< UnitCellmake_lattice_points (const Lattice &tiling_unit, const Lattice &superlattice, double tol)
 
double volume (const Lattice &lat)
 Returns the volume of a Lattice. More...
 
std::pair< bool, Eigen::Matrix3d > is_superlattice (const Lattice &scel, const Lattice &unit, double tol)
 Check if scel is a superlattice of unitcell unit and some integer transformation matrix T. More...
 
std::istream & operator>> (std::istream &in, const Lattice &lattice_in)
 
Lattice make_superduperlattice (const Lattice &lat1, const Lattice &lat2)
 returns Lattice that is smallest possible superlattice of both input Lattice More...
 
Eigen::Matrix3d cart2frac (const Eigen::Ref< const Eigen::Matrix3d > &cart_mat, const Lattice &lat)
 Returns 'frac_mat' which is transformation of 'cart_mat' if cart_vec_after = cart_mat*cart_vec then frac_vec_after = frac_mat*frac_vec where cart_vec = lat.lat_column_mat()*frac_vec and cart_vec_after = lat.lat_column_mat()*frac_vec_after. More...
 
Eigen::Matrix3d frac2cart (const Eigen::Ref< const Eigen::Matrix3d > &frac_mat, const Lattice &lat)
 Returns 'cart_mat' which is transformation of 'frac_mat' if cart_vec_after = cart_mat*cart_vec then frac_vec_after = frac_mat*frac_vec where cart_vec = lat.lat_column_mat()*frac_vec and cart_vec_after = lat.lat_column_mat()*frac_vec_after. More...
 
template<typename IntegralType , int Options = 0>
Lattice make_superlattice (const Lattice &lat, const Eigen::Matrix< IntegralType, 3, 3, Options > &transf_mat)
 Returns a super Lattice. Transformation matrix must be integer. More...
 
Eigen::Matrix3l make_transformation_matrix_to_super (const Lattice &tiling_unit, const Lattice &superlattice, double tol)
 
bool is_equivalent (const Lattice &ref_lattice, const Lattice &other)
 Check if ref_lattice = other*U, where U is unimodular. More...
 
bool compare_type (AtomPosition const &A, AtomPosition const &B, double tol)
 
bool operator== (Molecule const &A, Molecule const &B)
 
bool is_vacancy (const std::string &name)
 A vacancy is any Specie/Molecule with (name == "VA" || name == "va" || name == "Va") More...
 
bool is_molecule_name (const Molecule &mol, std::string name)
 Return true if Molecule name matches 'name', including Va checks. More...
 
bool less (LatticeNode const &A, LatticeNode const &B, double cost_tol)
 Compare two LatticeMap objects, based on their mapping cost first, followed by PrimGrid transformation matrices. More...
 
bool identical (LatticeNode const &A, LatticeNode const &B, double cost_tol)
 returns true if cost values and parent/child supercell transformations are same for A and B More...
 
bool identical (AssignmentNode const &A, AssignmentNode const &B)
 true if time_reversal and translation are identical More...
 
Eigen::Matrix3d const & get_matrix (MappingNode const &_node)
 External accessor for isometry, to provide xtal::SymOp adaptability. More...
 
Eigen::Vector3d const & get_translation (MappingNode const &_node)
 External accessor for translation, to provide xtal::SymOp adaptability. More...
 
bool get_time_reversal (MappingNode const &_node)
 External accessor for time_reversal, to provide xtal::SymOp adaptability. More...
 
Eigen::Matrix3i enforce_min_volume (const Lattice &unit, const Eigen::Matrix3i &T, const SymOpVector &point_grp, Index volume, bool fix_shape=false)
 Return a transformation matrix that ensures a supercell of at least some volume. More...
 
template<typename ExternSymGroupTypeIt >
Eigen::Matrix3i enforce_min_volume (ExternSymGroupTypeIt begin, ExternSymGroupTypeIt end, const Lattice &unit, const Eigen::Matrix3i &T, Index volume, bool fix_shape=false)
 
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...
 
template<typename ExternSymGroupTypeIt >
Eigen::Matrix3i canonical_hnf (ExternSymGroupTypeIt begin, ExternSymGroupTypeIt end, const Eigen::Matrix3i &T, const Lattice &ref_lattice)
 
std::vector< Indexinvariant_subgroup_indices (const Lattice &lat, SymOpVector const &super_grp)
 Construct indices of the subgroup that leaves a lattice unchanged. More...
 
template<typename ExternSymOpVector >
std::vector< Indexinvariant_subgroup_indices (const Lattice &lat, ExternSymOpVector const &super_grp)
 
template<typename OutputIt >
OutputIt invariant_subgroup_indices (const Lattice &lat, const std::vector< SymOp > &super_group, OutputIt result)
 Construct indices of the subgroup for which this->is_equivalent(copy_apply(op, *this)) More...
 
Lattice symmetrize (const Lattice &lat, const std::vector< SymOp > &enforced_group)
 Convert a Cartesian symmetry operation representation to fractional. More...
 
Lattice symmetrize (const Lattice &lat, double point_group_tolerance)
 Return a copy of the given lattice, which obeys the symmetry of its point group, when generated within the tolerance. More...
 
std::vector< SymOpmake_point_group (Lattice const &_lat)
 Populate. More...
 
std::vector< SymOpmake_point_group (Lattice const &_lat, double _tol)
 Populate. More...
 
template<typename Object , typename OpIterator >
std::pair< OpIterator, Eigen::Matrix3d > is_equivalent_superlattice (const Object &scel, const Object &unit, OpIterator begin, OpIterator end, double tol)
 
template<typename LatIterator , typename SymOpIterator >
Lattice make_equivalent_superduperlattice (LatIterator begin, LatIterator end, SymOpIterator op_begin, SymOpIterator op_end)
 returns Lattice that is smallest possible superlattice of all input Lattice More...
 
template<typename StructureType , typename ExternSymOpVector >
StructureType symmetrize (const StructureType &structure, const ExternSymOpVector &enforced_group)
 
SymOp operator* (const SymOp &LHS, const SymOp &RHS)
 Get a new SymOp that is equivalent to subsequent application of both SymOps. More...
 
const SymOpMatrixTypeget_matrix (const SymOp &op)
 Accessor for SymOpType. Returns transformation matrix (Cartesian). More...
 
const SymOpTranslationTypeget_translation (const SymOp &op)
 Accessor for SymOpType. Returns translation vector (tau). More...
 
SymOpTimeReversalType get_time_reversal (const SymOp &op)
 
template<typename SymOpCompareType >
void close_group (SymOpVector *partial_group, const SymOpCompareType &symop_binary_comp)
 
template<typename SymOpCompareType , typename... CompareArgs>
void close_group (SymOpVector *partial_group, const CompareArgs &... args)
 
Coordinate make_superlattice_coordinate (const UnitCell &ijk, const Superlattice &superlattice)
 
Coordinate make_superlattice_coordinate (const UnitCell &ijk, const Lattice &tiling_unit, const Lattice &superlattice)
 
template BasicStructure make_superstructure< int, 0 > (const BasicStructure &tiling_unit, const Eigen::Matrix< int, 3, 3, 0 > &transformation_matrix)
 
template BasicStructure make_superstructure< long, 0 > (const BasicStructure &tiling_unit, const Eigen::Matrix< long, 3, 3, 0 > &transformation_matrix)
 
std::set< Eigen::Matrix3d, StandardOrientationCompare_niggli_set (const Lattice &in_lat, double compare_tol, bool keep_handedness)
 
void _formatted_print (std::ostream &stream, Eigen::Vector3d vec, COORD_TYPE mode, char term, Eigen::IOFormat format)
 
std::ostream & operator<< (std::ostream &stream, const Coordinate &coord)
 
std::vector< Eigen::Matrix3d > _skew_transforms ()
 
Lattice replace_vector (const Lattice &lat, const Eigen::Vector3d &new_vector, double tol)
 
std::vector< Eigen::Matrix3d > _cell_invariant_transforms ()
 
bool standard_orientation_spatial_compare (const Eigen::Matrix3d &low_score_lat_mat, const Eigen::Matrix3d &high_score_lat_mat, double compare_tol)
 
std::vector< Indexinvariant_subgroup_indices (const Lattice &lat, std::vector< SymOp > const &super_grp)
 Construct the subgroup that leaves a lattice unchanged. More...
 
Index niggli_index (const Eigen::Matrix3d &test_lat_mat, double compare_tol)
 Number of niggli criteria met. More...
 
Lattice niggli (const Lattice &in_lat, double compare_tol, bool keep_handedness=false)
 
bool is_niggli (const Eigen::Matrix3d &test_lat_mat, double compare_tol)
 
bool is_niggli (const Lattice &test_lat, double compare_tol)
 Check whether the given lattice is niggli (does not check for orientation) More...
 
Eigen::VectorXd spatial_unroll (const Eigen::Matrix3d &lat_mat, double compare_tol)
 Generate a vector whose lexicographical value determines how well it's oriented in space. More...
 
bool standard_orientation_spatial_compare (const Eigen::Matrix3d &low_score_lat_mat, Eigen::Matrix3d &high_score_lat_mat, double compare_tol)
 Compare the spatial orientation (ignoring matrix symmetry) and determine which one is oriented more standard. True if high is more standard. More...
 
bool standard_orientation_compare (const Eigen::Matrix3d &low_score_lat_mat, const Eigen::Matrix3d &high_score_lat_mat, double compare_tol)
 Determine whether high_score has a more standard format than low_score. More...
 
SimpleStructure make_superstructure (Eigen::Ref< const Eigen::Matrix3i > const &_T, SimpleStructure const &_sstruc)
 
std::vector< Indexsuperstructure_basis_idx (Eigen::Ref< const Eigen::Matrix3i > const &_T, SimpleStructure const &_sstruc)
 Constructs a vector containing the basis index of the ith site in the supercell. More...
 
SimpleStructure make_simple_structure (BasicStructure const &_struc)
 Construct from decorated structure. More...
 
std::vector< std::set< Index > > atom_site_compatibility (SimpleStructure const &sstruc, BasicStructure const &_prim)
 Determine which sites of a BasicStructure can host each atom of a SimpleStructure result[i] is set of site indices in. More...
 
std::vector< std::set< Index > > mol_site_compatibility (SimpleStructure const &sstruc, BasicStructure const &_prim)
 Determine which sites of a BasicStructure can host each molecule of a SimpleStructure result[i] is set of site indices in. More...
 
void _atomize (SimpleStructure &_sstruc, Eigen::Ref< const Eigen::VectorXi > const &_mol_occ, BasicStructure const &_reference)
 use information in _reference to initialize atom_info from mol_info More...
 
BasicStructure make_basic_structure (SimpleStructure const &_sstruc, std::vector< DoFKey > const &_all_dofs, SimpleStructureTools::SpeciesMode mode, std::vector< std::vector< Molecule >> _allowed_occupants={})
 Construct BasicStructure from SimpleStructure. More...
 
std::vector< Eigen::MatrixXd > generate_invariant_shuffle_modes (const std::vector< xtal::SymOp > &factor_group, const std::vector< Eigen::PermutationMatrix< Eigen::Dynamic, Eigen::Dynamic, Index >> &permute_group)
 
std::ostream & operator<< (std::ostream &stream, const Site &site)
 
Site operator* (const SymOp &LHS, const Site &RHS)
 Apply the symmetry operation to the site, and return a new transformed value. More...
 
template<typename ExternSymOp >
Site operator* (const ExternSymOp &LHS, const Site &RHS)
 
Site operator+ (const Site &LHS, const Coordinate &RHS)
 
Site operator+ (const Coordinate &LHS, const Site &RHS)
 
std::ostream & operator<< (std::ostream &sout, const UnitCellCoord &site)
 

Typedef Documentation

◆ LinearIndexConverter

template<typename IntegralCoordinateType >
using CASM::xtal::LinearIndexConverter = typedef typename _LinearIndexConverter<IntegralCoordinateType>::type

Definition at line 186 of file LinearIndexConverter.hh.

◆ SymOpMatrixType

typedef Eigen::Matrix3d CASM::xtal::SymOpMatrixType

Definition at line 15 of file SymType.hh.

◆ SymOpTimeReversalType

Definition at line 17 of file SymType.hh.

◆ SymOpTranslationType

typedef Eigen::Vector3d CASM::xtal::SymOpTranslationType

Definition at line 16 of file SymType.hh.

◆ SymOpVector

typedef std::vector< SymOp > CASM::xtal::SymOpVector

This defines the type of the object representing symmetry operations within the crystallography classes. Any symmetry related operations within the crystallography module must be in terms of this type.

Definition at line 13 of file CanonicalForm.hh.

Function Documentation

◆ _atomize()

void CASM::xtal::_atomize ( SimpleStructure _sstruc,
Eigen::Ref< const Eigen::VectorXi > const &  _mol_occ,
BasicStructure const &  _reference 
)

use information in _reference to initialize atom_info from mol_info

Definition at line 258 of file SimpleStructureTools.cc.

◆ _cell_invariant_transforms()

std::vector<Eigen::Matrix3d> CASM::xtal::_cell_invariant_transforms ( )

Definition at line 12 of file Niggli.cc.

◆ _formatted_print()

void CASM::xtal::_formatted_print ( std::ostream &  stream,
Eigen::Vector3d  vec,
COORD_TYPE  mode,
char  term,
Eigen::IOFormat  format 
)

Definition at line 109 of file Coordinate.cc.

◆ _niggli_set()

std::set< Eigen::Matrix3d, StandardOrientationCompare > CASM::xtal::_niggli_set ( const Lattice in_lat,
double  compare_tol,
bool  keep_handedness 
)

Definition at line 284 of file Niggli.cc.

◆ _skew_transforms()

std::vector<Eigen::Matrix3d> CASM::xtal::_skew_transforms ( )

Definition at line 26 of file Lattice.cc.

◆ atom_site_compatibility()

std::vector< std::set< Index > > CASM::xtal::atom_site_compatibility ( SimpleStructure const &  sstruc,
BasicStructure const &  _prim 
)

Determine which sites of a BasicStructure can host each atom of a SimpleStructure result[i] is set of site indices in.

Parameters
_primthat can host atom 'i' of
sstruc

Definition at line 348 of file SimpleStructureTools.cc.

◆ calc_rotation_angle_and_axis()

std::pair< double, Eigen::Vector3d > CASM::xtal::calc_rotation_angle_and_axis ( const SymOp op,
const Lattice lat 
)

Calculates the rotation angle and axis of a symmetry operation. This function is almost exactly identical to the constructor of SymInfo::SymInfo.

Calculates the rotation angle and axis of a symmetry operation. This function is almost exactly identical to the constructor of SymInfo::SymInfo

Definition at line 323 of file BasicStructureTools.cc.

◆ close_group() [1/2]

template<typename SymOpCompareType , typename... CompareArgs>
void CASM::xtal::close_group ( SymOpVector partial_group,
const CompareArgs &...  args 
)

Definition at line 128 of file SymType.hh.

◆ close_group() [2/2]

template<typename SymOpCompareType >
void CASM::xtal::close_group ( SymOpVector partial_group,
const SymOpCompareType &  symop_binary_comp 
)

Combines every pair of symmetry operations and adds any missing resulting operations to the group. In order to determine what a missing operation is, a comparator type must be provided, and any arguments needed to construct it must be passed as final arguments.

For example, when closing a factor group, we consider operations to be the same when the transformation matrix matches to within a tolerance, and the translation components are equivalent under translational symmetry. We therefore want to use the SymOpPeriodicCompare_f comparator, which requires a lattice and a tolerance. The group closure call for this situation might look like:

close_group<SymOpPeriodicCompare_f> close_group(&my_partial_group, my_lattice, CASM::TOL);

Definition at line 89 of file SymType.hh.

◆ component_descriptions()

std::vector< std::string > CASM::xtal::component_descriptions ( DoFSet const &  dofset)

Returns descriptive names of the components in a DoFSet, using AnisoValTraits::variable_descriptors()

Definition at line 12 of file DoFSet.cc.

◆ dofset_transformation_matrix()

Eigen::MatrixXd CASM::xtal::dofset_transformation_matrix ( const Eigen::MatrixXd &  from_basis,
const Eigen::MatrixXd &  to_basis,
double  tol 
)

Create the symmtery representation for going from one basis to another.

Definition at line 52 of file DoFSet.cc.

◆ find_index()

Index CASM::xtal::find_index ( const std::vector< Site > &  basis,
const Site test_site,
double  tol 
)

return basis index of site that is same type and within distance 'tol' (in Angstr) of test_site If no such site exists in basis, return the size of the basis

Definition at line 233 of file BasicStructureTools.cc.

◆ generate_invariant_shuffle_modes()

std::vector<Eigen::MatrixXd> CASM::xtal::generate_invariant_shuffle_modes ( const std::vector< xtal::SymOp > &  factor_group,
const std::vector< Eigen::PermutationMatrix< Eigen::Dynamic, Eigen::Dynamic, Index >> &  permute_group 
)

Calculates the invariant shuffle modes in the primitive unit cell. They symmetry preserving distortions are found by applying the Reynolds operator to the basis of displacement vectors. The average of the resulting basis vectors is used to form an orthonormal basis.

Parameters
factor_groupuse make_factor_group(struc) to obtain this group
permute_groupuse make_permute_group(struc,factor_group) to obtain this group
Returns
the vector of shuffle modes that are invariant under symmetry. Each element has a size 3 x basis_size.

Definition at line 191 of file SimpleStructureTools.cc.

◆ get_matrix() [1/2]

const SymOpMatrixType & CASM::xtal::get_matrix ( const SymOp op)

Accessor for SymOpType. Returns transformation matrix (Cartesian).

Definition at line 17 of file SymType.cc.

◆ get_matrix() [2/2]

Eigen::Matrix3d const& CASM::xtal::get_matrix ( MappingNode const &  _node)
inline

External accessor for isometry, to provide xtal::SymOp adaptability.

Definition at line 382 of file StrucMapping.hh.

◆ get_time_reversal() [1/2]

SymOpTimeReversalType CASM::xtal::get_time_reversal ( const SymOp op)

Accessor for SymOpType. Returns whether the symmetry operation is time reversal active.

Definition at line 23 of file SymType.cc.

◆ get_time_reversal() [2/2]

bool CASM::xtal::get_time_reversal ( MappingNode const &  _node)
inline

External accessor for time_reversal, to provide xtal::SymOp adaptability.

Definition at line 394 of file StrucMapping.hh.

◆ get_translation() [1/2]

const SymOpTranslationType & CASM::xtal::get_translation ( const SymOp op)

Accessor for SymOpType. Returns translation vector (tau).

Definition at line 19 of file SymType.cc.

◆ get_translation() [2/2]

Eigen::Vector3d const& CASM::xtal::get_translation ( MappingNode const &  _node)
inline

External accessor for translation, to provide xtal::SymOp adaptability.

Definition at line 388 of file StrucMapping.hh.

◆ identical() [1/2]

bool CASM::xtal::identical ( AssignmentNode const &  A,
AssignmentNode const &  B 
)

true if time_reversal and translation are identical

Definition at line 333 of file StrucMapping.cc.

◆ identical() [2/2]

bool CASM::xtal::identical ( LatticeNode const &  A,
LatticeNode const &  B,
double  cost_tol 
)

returns true if cost values and parent/child supercell transformations are same for A and B

Definition at line 310 of file StrucMapping.cc.

◆ invariant_subgroup_indices() [1/4]

template<typename OutputIt >
OutputIt CASM::xtal::invariant_subgroup_indices ( const Lattice lat,
const std::vector< SymOp > &  super_group,
OutputIt  result 
)

Construct indices of the subgroup for which this->is_equivalent(copy_apply(op, *this))

Definition at line 45 of file SymTools.hh.

◆ invariant_subgroup_indices() [2/4]

template<typename ExternSymOpVector >
std::vector<Index> CASM::xtal::invariant_subgroup_indices ( const Lattice lat,
ExternSymOpVector const &  super_grp 
)

Definition at line 36 of file SymTools.hh.

◆ invariant_subgroup_indices() [3/4]

std::vector<Index> CASM::xtal::invariant_subgroup_indices ( const Lattice lat,
std::vector< SymOp > const &  super_grp 
)

Construct the subgroup that leaves a lattice unchanged.

Definition at line 88 of file SymTools.cc.

◆ invariant_subgroup_indices() [4/4]

std::vector<Index> CASM::xtal::invariant_subgroup_indices ( const Lattice lat,
SymOpVector const &  super_grp 
)

Construct indices of the subgroup that leaves a lattice unchanged.

◆ is_equivalent()

bool CASM::xtal::is_equivalent ( const Lattice ref_lattice,
const Lattice other 
)

Check if ref_lattice = other*U, where U is unimodular.

Definition at line 100 of file LatticeIsEquivalent.cc.

◆ is_equivalent_superlattice()

template<typename Object , typename OpIterator >
std::pair<OpIterator, Eigen::Matrix3d> CASM::xtal::is_equivalent_superlattice ( const Object &  scel,
const Object &  unit,
OpIterator  begin,
OpIterator  end,
double  tol 
)

Check if there is a symmetry operation, op, and transformation matrix T, such that scel is a superlattice of the result of applying op to unit

Returns
pair corresponding to first successful op and T, or with op=end if not successful

Definition at line 110 of file SymTools.hh.

◆ is_niggli() [1/2]

bool CASM::xtal::is_niggli ( const Eigen::Matrix3d &  test_lat_mat,
double  compare_tol 
)

Check whether the given lattice (represented as a matrix) is in niggli TYPE ?? reduced form (does not check for orientation)

Definition at line 271 of file Niggli.cc.

◆ is_niggli() [2/2]

bool CASM::xtal::is_niggli ( const Lattice test_lat,
double  compare_tol 
)

Check whether the given lattice is niggli (does not check for orientation)

Definition at line 278 of file Niggli.cc.

◆ is_primitive()

bool CASM::xtal::is_primitive ( const BasicStructure struc,
double  tol = TOL 
)

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

Definition at line 244 of file BasicStructureTools.cc.

◆ less()

bool CASM::xtal::less ( LatticeNode const &  A,
LatticeNode const &  B,
double  cost_tol 
)

Compare two LatticeMap objects, based on their mapping cost first, followed by PrimGrid transformation matrices.

Compare two LatticeMap objects, based on their mapping cost.

Definition at line 295 of file StrucMapping.cc.

◆ make_asymmetric_unit() [1/2]

std::set< std::set< Index > > CASM::xtal::make_asymmetric_unit ( const xtal::BasicStructure struc)

Return indices of equivalent basis sites.

Return indices of equivalent basis sites

Equivalent to:

std::vector<SymOp> factor_group = make_factor_group(struc,
struc.lattice().tol());
SharedPrimFormatter< jsonParser > factor_group()
std::set< std::set< Index > > make_asymmetric_unit(const xtal::BasicStructure &struc, const std::vector< SymOp > &factor_group)
Return indices of equivalent basis sites.
std::vector< SymOp > make_factor_group(const BasicStructure &struc, double tol=TOL)

Definition at line 602 of file BasicStructureTools.cc.

◆ make_asymmetric_unit() [2/2]

std::set< std::set< Index > > CASM::xtal::make_asymmetric_unit ( const xtal::BasicStructure struc,
const std::vector< SymOp > &  factor_group 
)

Return indices of equivalent basis sites.

Return indices of equivalent basis sites

Parameters
strucInput structure
factor_groupSymmetry operations of the input structure.
Returns
One std::set<Index> for each orbit of equivalent basis sites.

Note:

  • This function assumes the operations in factor_group do in fact map sites onto equivalents sites and does not double-check site equivalence.

Definition at line 573 of file BasicStructureTools.cc.

◆ make_basic_structure()

BasicStructure CASM::xtal::make_basic_structure ( SimpleStructure const &  _sstruc,
std::vector< DoFKey > const &  _all_dofs,
SimpleStructureTools::SpeciesMode  mode,
std::vector< std::vector< Molecule >>  _allowed_occupants = {} 
)

Construct BasicStructure from SimpleStructure.

Parameters
_sstrucSimpleStructure used as source data for conversion
_all_dofsholds names of additional DOFs to initialize in structure
modespecifies whether ATOM or MOL info of _sstruc should be used to build sites of structure
_allowed_occupantsList of allowed molecules at each site; if empty, occupants are assumed to be atoms having the species names and attributes indicated by _sstruc

Definition at line 118 of file SimpleStructureTools.cc.

◆ make_bring_within_f()

IntegralCoordinateWithin_f CASM::xtal::make_bring_within_f ( const Supercell scel)

Make IntegralCoordinateWithin_f for Supercell [deprecated].

Definition at line 520 of file Supercell.cc.

◆ make_dofset_map()

template<typename DoFSetType >
std::map<std::string, DoFSetType> CASM::xtal::make_dofset_map ( std::vector< DoFSetType > const &  dofset_vec)

Definition at line 130 of file DoFSet.hh.

◆ make_equivalent_superduperlattice()

template<typename LatIterator , typename SymOpIterator >
Lattice CASM::xtal::make_equivalent_superduperlattice ( LatIterator  begin,
LatIterator  end,
SymOpIterator  op_begin,
SymOpIterator  op_end 
)

returns Lattice that is smallest possible superlattice of all input Lattice

SymOpIterator are provided to apply to each Lattice in an attempt to find the smallest possible superduperlattice of all symmetrically transformed Lattice

Definition at line 130 of file SymTools.hh.

◆ make_factor_group()

std::vector< SymOp > CASM::xtal::make_factor_group ( const BasicStructure struc,
double  tol = TOL 
)

Create the factor group of the given structure. If the structure has no degrees of freedom affected by time reversal, time reversal is ignored. Otherwise symmetry operations are checked for time reversal

Definition at line 466 of file BasicStructureTools.cc.

◆ make_lattice_points() [1/2]

std::vector< UnitCell > CASM::xtal::make_lattice_points ( const Eigen::Matrix3l transformation_matrix)

Definition at line 108 of file IntegralCoordinateWithin.cc.

◆ make_lattice_points() [2/2]

std::vector< UnitCell > CASM::xtal::make_lattice_points ( const Lattice tiling_unit,
const Lattice superlattice,
double  tol 
)

Returns all the lattice points that exists when tiling the tiling unit inside the superlattice

Definition at line 113 of file IntegralCoordinateWithin.cc.

◆ make_permutation_representation()

std::vector< Eigen::PermutationMatrix< Eigen::Dynamic, Eigen::Dynamic, Index > > CASM::xtal::make_permutation_representation ( const xtal::BasicStructure struc,
const std::vector< SymOp > &  factor_group 
)

Create the permutation group of a structure.

Create the permutation group of a structure.

The permutation group can be used to identify how the basis of a structure is transformed under the application of a symmetry operation. The permutation representation for a structure is obtained by applying each factor group operation to the structure and identifying the index of the symmetry transformed site.

Parameters
strucThe input structure
factor_groupThe factor group of struc, such as generated with xtal::make_factor_group(struc).
Returns
The permutation matrix at index idx corresponds to the effect of the symmetry operation in factor_group[idx] on the basis.

The Eigen::PermutationMatrix value is stored in a vector of indices():

  • values_after = values_before * permutation_matrix;
  • values_after[permutation_matrix.indices()[i]] = values_before[i]

Definition at line 523 of file BasicStructureTools.cc.

◆ make_point_group() [1/2]

std::vector< SymOp > CASM::xtal::make_point_group ( Lattice const &  _lat)

Populate.

Relax the vectors of the given lattice such that it obeys the symmetry of the given group, where the symmetry operations are given in fractional representations

Parameters
point_groupwith the point group of this lattice
point_groupshould be empty
pg_tolcan be increased to find point group of lattice vectors that are slightly distorted due to numerical noise

Definition at line 113 of file SymTools.cc.

◆ make_point_group() [2/2]

std::vector< SymOp > CASM::xtal::make_point_group ( Lattice const &  _lat,
double  _tol 
)

Populate.

Parameters
point_groupwith the point group of this lattice
point_groupshould be empty
pg_tolcan be increased to find point group of lattice vectors that are slightly distorted due to numerical noise

Definition at line 117 of file SymTools.cc.

◆ make_primitive()

BasicStructure CASM::xtal::make_primitive ( const BasicStructure non_primitive_struc,
double  tol = TOL 
)

Returns the smallest possible tiling unit of the given structure.

Definition at line 250 of file BasicStructureTools.cc.

◆ make_simple_structure()

SimpleStructure CASM::xtal::make_simple_structure ( BasicStructure const &  _struc)

Construct from decorated structure.

Definition at line 99 of file SimpleStructureTools.cc.

◆ make_superlattice_coordinate() [1/2]

Coordinate CASM::xtal::make_superlattice_coordinate ( const UnitCell ijk,
const Lattice tiling_unit,
const Lattice superlattice 
)

Definition at line 104 of file UnitCellCoord.cc.

◆ make_superlattice_coordinate() [2/2]

Coordinate CASM::xtal::make_superlattice_coordinate ( const UnitCell ijk,
const Superlattice superlattice 
)

Converts the position of a lattice site into a Coordinate within the superlattice

Definition at line 112 of file UnitCellCoord.cc.

◆ make_superstructure() [1/2]

template<typename IntegralType , int Options = 0>
BasicStructure CASM::xtal::make_superstructure ( const BasicStructure tiling_unit,
const Eigen::Matrix< IntegralType, 3, 3, Options > &  transformation_matrix 
)

Given an integer transformation matrix, create a superstructure whose lattice is the product of the original lattice and the transformation matrix.

Definition at line 653 of file BasicStructureTools.cc.

◆ make_superstructure() [2/2]

SimpleStructure CASM::xtal::make_superstructure ( Eigen::Ref< const Eigen::Matrix3i > const &  _T,
SimpleStructure const &  _sstruc 
)

Definition at line 45 of file SimpleStructureTools.cc.

◆ make_superstructure< int, 0 >()

template BasicStructure CASM::xtal::make_superstructure< int, 0 > ( const BasicStructure tiling_unit,
const Eigen::Matrix< int, 3, 3, 0 > &  transformation_matrix 
)

◆ make_superstructure< long, 0 >()

template BasicStructure CASM::xtal::make_superstructure< long, 0 > ( const BasicStructure tiling_unit,
const Eigen::Matrix< long, 3, 3, 0 > &  transformation_matrix 
)

◆ make_transformation_matrix_to_super()

Eigen::Matrix3l CASM::xtal::make_transformation_matrix_to_super ( const Lattice tiling_unit,
const Lattice superlattice,
double  tol 
)

Calculates the transformation matrix that takes the tiling unit to the superlattice. Throws exceptions if the superlattice isn't compatible with its tiling unit

Definition at line 848 of file Lattice.cc.

◆ mol_site_compatibility()

std::vector< std::set< Index > > CASM::xtal::mol_site_compatibility ( SimpleStructure const &  sstruc,
BasicStructure const &  _prim 
)

Determine which sites of a BasicStructure can host each molecule of a SimpleStructure result[i] is set of site indices in.

Parameters
_primthat can host molecule 'i' of
sstruc

Definition at line 333 of file SimpleStructureTools.cc.

◆ niggli()

Lattice CASM::xtal::niggli ( const Lattice in_lat,
double  compare_tol,
bool  keep_handedness = false 
)

Convert the given lattice into it's niggli reduced form, with the most standard orientation possilbe

Converts a given Lattice into a new Lattice that satisfies the Niggli criteria. Because there are several representations of the Lattice that satisfy these conditions, the one with the most standard orientation will be returned.

The returned Niggli cell will be right handed regardless of the inputted lattice, unless keep_handedness is set to true.

Definition at line 330 of file Niggli.cc.

◆ niggli_index()

Index CASM::xtal::niggli_index ( const Eigen::Matrix3d &  test_lat_mat,
double  compare_tol 
)

Number of niggli criteria met.

Definition at line 267 of file Niggli.cc.

◆ operator*() [1/3]

template<typename ExternSymOp >
Site CASM::xtal::operator* ( const ExternSymOp &  LHS,
const Site RHS 
)

Definition at line 128 of file Site.hh.

◆ operator*() [2/3]

Site CASM::xtal::operator* ( const SymOp LHS,
const Site RHS 
)

Apply the symmetry operation to the site, and return a new transformed value.

Definition at line 456 of file Site.cc.

◆ operator*() [3/3]

SymOp CASM::xtal::operator* ( const SymOp LHS,
const SymOp RHS 
)

Get a new SymOp that is equivalent to subsequent application of both SymOps.

Definition at line 8 of file SymType.cc.

◆ operator+() [1/2]

Site CASM::xtal::operator+ ( const Coordinate LHS,
const Site RHS 
)

Definition at line 468 of file Site.cc.

◆ operator+() [2/2]

Site CASM::xtal::operator+ ( const Site LHS,
const Coordinate RHS 
)

Definition at line 462 of file Site.cc.

◆ operator<<() [1/3]

std::ostream& CASM::xtal::operator<< ( std::ostream &  sout,
const UnitCellCoord site 
)
inline

Definition at line 191 of file UnitCellCoord.hh.

◆ operator<<() [2/3]

std::ostream& CASM::xtal::operator<< ( std::ostream &  stream,
const Coordinate coord 
)

Definition at line 339 of file Coordinate.cc.

◆ operator<<() [3/3]

std::ostream & CASM::xtal::operator<< ( std::ostream &  stream,
const Site site 
)

Definition at line 417 of file Site.cc.

◆ replace_vector()

Lattice CASM::xtal::replace_vector ( const Lattice lat,
const Eigen::Vector3d &  new_vector,
double  tol 
)
related

Definition at line 812 of file Lattice.cc.

◆ sort_factor_group()

void CASM::xtal::sort_factor_group ( std::vector< SymOp > &  factor_group,
const Lattice lat 
)

Sort SymOp in the SymGroup.

Sort a factor group based on a lexicographical comparison of (-det, -trace, angle, axis, tau)

SymOp are sorted by lexicographical comparison of: (-det, -trace, angle, axis, tau)

  • angle is positive
  • axis[0] is positive

This function is an almost identical clone of SymGroup::sort()

Definition at line 383 of file BasicStructureTools.cc.

◆ spatial_unroll()

Eigen::VectorXd CASM::xtal::spatial_unroll ( const Eigen::Matrix3d &  lat_mat,
double  compare_tol 
)

Generate a vector whose lexicographical value determines how well it's oriented in space.

Generates a vector of values from a Lattice column matrix, whose lexicographical value determines how standard it's orientation is.

First entry is the number of non-negative entries in the matrix. Next entries are the diagonal values. Next entries are the negative of the absolute value of the off diagonal entries. Finally, the signs of the off diagonal entries.

Maximizing the lexicographical value will favor matrices that are diagonal if possible, followed by a preference for upper triangular matrices. Off diagonal entries are minimized.

Definition at line 354 of file Niggli.cc.

◆ standard_orientation_compare()

bool CASM::xtal::standard_orientation_compare ( const Eigen::Matrix3d &  low_score_lat_mat,
const Eigen::Matrix3d &  high_score_lat_mat,
double  compare_tol 
)

Determine whether high_score has a more standard format than low_score.

Compares two lattice column matrices to each other and determines which one has a more standard orientation based on whether said matrices are (bi)symmetric and how aligned the vectors are along the Cartesian directions.

A bisymmetric matrix will always have better standard orientation than one that is simply symmetric. In turn, a symmetric matrix will always have a better standard orientation than on that is not. The criteria with lowest impact is the orientation of the lattice vectors in space, as defined in standard_orientation_spatial_compare

The routine returns "low_score_mat<high_score_mat", which is true if high_score_mat has a more standard orientation.

Definition at line 427 of file Niggli.cc.

◆ standard_orientation_spatial_compare() [1/2]

bool CASM::xtal::standard_orientation_spatial_compare ( const Eigen::Matrix3d &  low_score_lat_mat,
const Eigen::Matrix3d &  high_score_lat_mat,
double  compare_tol 
)

Given two lattice column matrices, determine which one of them is oriented in a more standard manner in space relative to the Cartesian coordinates.

The routine returns "low_score_lat_mat<high_score_lat_mat", which is true if high_score_lat_mat has a more standard spatial orientation.

Definition at line 400 of file Niggli.cc.

◆ standard_orientation_spatial_compare() [2/2]

bool CASM::xtal::standard_orientation_spatial_compare ( const Eigen::Matrix3d &  low_score_lat_mat,
Eigen::Matrix3d &  high_score_lat_mat,
double  compare_tol 
)

Compare the spatial orientation (ignoring matrix symmetry) and determine which one is oriented more standard. True if high is more standard.

◆ superstructure_basis_idx()

std::vector< Index > CASM::xtal::superstructure_basis_idx ( Eigen::Ref< const Eigen::Matrix3i > const &  _T,
SimpleStructure const &  _sstruc 
)

Constructs a vector containing the basis index of the ith site in the supercell.

Calculates the parent basis index of each site in a supercell that is generated with the make_superstructure method

Parameters
_Tis the transformation matrix linking _sstruc and the supercell
_sstrucis a SimpleStructure
Returns
std::vector<Index> with the Index in the ith entry corresponding to the index of site i in _sstruc

Definition at line 84 of file SimpleStructureTools.cc.

◆ symmetrize() [1/4]

BasicStructure CASM::xtal::symmetrize ( const BasicStructure structure,
const std::vector< SymOp > &  enforced_group 
)

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.

Definition at line 609 of file BasicStructureTools.cc.

◆ symmetrize() [2/4]

Lattice CASM::xtal::symmetrize ( const Lattice lat,
const std::vector< SymOp > &  enforced_group 
)

Convert a Cartesian symmetry operation representation to fractional.

Convert a fractional symmetry operation representation to Cartesian

Return a copy of the given lattice, which obeys the symmetry of the given group

Parameters
enforced_group

Definition at line 95 of file SymTools.cc.

◆ symmetrize() [3/4]

Lattice CASM::xtal::symmetrize ( const Lattice lat,
double  point_group_tolerance 
)

Return a copy of the given lattice, which obeys the symmetry of its point group, when generated within the tolerance.

Parameters
point_group_tolerance

Definition at line 108 of file SymTools.cc.

◆ symmetrize() [4/4]

template<typename StructureType , typename ExternSymOpVector >
StructureType CASM::xtal::symmetrize ( const StructureType &  structure,
const ExternSymOpVector &  enforced_group 
)

Definition at line 155 of file SymTools.hh.