CASM
1.1.0
A Clusters Approach to Statistical Mechanics
|
#include <ConfigDoF.hh>
A container class for the different degrees of freedom values a Configuration might have
ConfigDoF has:
The values of the discrete site DoF within the supercell ("occupation", Eigen::VectorXi)
Example: Occupation values, accessed via Eigen::VectorXi const &ConfigDoF::occupation()
, with integer value corresponding to which Molecule in the the a Site::occupant_dof vector is occupying a particular site:
[<- sublattice 0 "occ" values -> | <- sublattice 1 "occ" values -> | ... ]
The values of the continuous site DoF within the supercell ("local_dofs", std::map<DoFKey, LocalContinuousConfigDoFValues>)
Example: Displacement values, with prim DoF basis equal to the standard basis (dx, dy, dz), accessed via Eigen::MatrixXd const &ConfigDoF::local_dofs("disp").values()
:
[<- sublattice 0 dx values -> | <- sublattice 1 dx values -> | ... ] [<- sublattice 0 dy values -> | <- sublattice 1 dy values -> | ... ] [<- sublattice 0 dz values -> | <- sublattice 1 dz values -> | ... ]
Example: Displacement values, with non-standard prim DoF basis:
"basis" : [ { "coordinate": [c0x, c0y, c0z], "occupants": [...], "dofs": { "disp" : { "axis_names" : ["dxy", "dz"], "axes" : [[1.0, 1.0, 0.0], [0.0, 0.0, 1.0]]}} }, { "coordinate": [c1x, c1y, c1z], "occupants": [...], "dofs": { "disp" : { "axis_names" : ["d\bar{x}y", "dz"], "axes" : [[-1.0, 1.0, 0.0], [0.0, 0.0, 1.0]]}}
... }
[<- sublat 0 dxy values -> | <- sublat 1 d\bar{x}y values ->| ... ] [<- sublat 0 dz values -> | <- sublat 1 dz values -> | ... ]
Note that the values matrix has only two rows, this is the maximum site basis dimension.
Example: Displacement values, with varying prim DoF basis:
"basis" : [ { "coordinate": [c0x, c0y, c0z], "occupants": [...], "dofs": { "disp" : {} }, { "coordinate": [c1x, c1y, c1z], "occupants": [...], "dofs": { "disp" : { "axis_names" : ["dxy", "dz"], "axes" : [[-1.0, 1.0, 0.0], [0.0, 0.0, 1.0]]}}
... }
[<- sublattice 0 dx values -> | <- sublattice 1 dxy values ->| ... ] [<- sublattice 0 dy values -> | <- sublattice 1 dz values -> | ... ] [<- sublattice 0 dz values -> | <- 0.0 values -> | ... ]
Note that the values matrix has three rows, this is the maximum site basis dimension, but for sublattices with lower site basis dimension it is padded with fixed zeros.
The values of the continuous global DoF ("global_dofs", std::map<DoFKey, GlobalContinuousConfigDoFValues>).
Example: GLstrain values, with prim DoF basis equal to the standard basis, accessed via Eigen::VectorXd const &ConfigDoF::global_dofs("GLstrain").values()
:
[e_xx, e_yy, e_zz, sqrt(2)*e_yz, sqrt(2)*e_xz, sqrt(2)*e_xy]
Note: Continuous DoF values stored in memory in ConfigDoF are coordinates in the prim DoF basis (the "axes" given in "prim.json" which set xtal::SiteDoFSet::basis() or xtal::DoFSet::basis()), but when saved to file (i.e. .casm/config/config_list.json
) they are saved as coordinates in the standard DoF basis (with axes meaning as described by AnisoValTraits::standard_var_names()).
Definition at line 117 of file ConfigDoF.hh.
Public Member Functions | |
template<typename GlobalInfoContainerType , typename LocalInfoContainerType > | |
ConfigDoF (Index _N_sublat, Index _N_vol, GlobalInfoContainerType const &global_dof_info, LocalInfoContainerType const &local_dof_info, std::vector< SymGroupRepID > const &occ_symrep_IDs, double _tol) | |
Index | size () const |
Number of sites in the ConfigDoF. More... | |
Index | n_vol () const |
Integer volume of ConfigDoF. More... | |
Index | n_sublat () const |
Number of sublattices in ConfigDoF. More... | |
double | tol () const |
Tolerance for comparison of continuous DoF values. More... | |
void | setZero () |
Set all DoF values to zero. More... | |
int & | occ (Index i) |
Reference occupation value on site i. More... | |
const int & | occ (Index i) const |
Const reference to occupation value on site i. More... | |
void | set_occupation (Eigen::Ref< const Eigen::VectorXi > const &_occupation) |
Set occupation values. More... | |
Eigen::VectorXi const & | occupation () const |
Const reference occupation values. More... | |
bool | has_occupation () const |
std::map< DoFKey, GlobalContinuousConfigDoFValues > const & | global_dofs () const |
GlobalContinuousConfigDoFValues const & | global_dof (DoFKey const &_key) const |
GlobalContinuousConfigDoFValues & | global_dof (DoFKey const &_key) |
bool | has_global_dof (DoFKey const &_key) const |
void | set_global_dof (DoFKey const &_key, Eigen::Ref< const Eigen::VectorXd > const &_val) |
Set global continuous DoF values. More... | |
std::map< DoFKey, LocalContinuousConfigDoFValues > const & | local_dofs () const |
LocalContinuousConfigDoFValues const & | local_dof (DoFKey const &_key) const |
LocalContinuousConfigDoFValues & | local_dof (DoFKey const &_key) |
bool | has_local_dof (DoFKey const &_key) const |
void | set_local_dof (DoFKey const &_key, Eigen::Ref< const Eigen::MatrixXd > const &_val) |
Set local continuous DoF values. More... | |
ConfigDoF & | apply_sym (PermuteIterator const &it) |
ConfigDoF & | apply_sym_no_permute (SymOp const &_op) |
void | swap (ConfigDoF &RHS) |
Private Attributes | |
LocalDiscreteConfigDoFValues | m_occupation |
std::map< std::string, GlobalContinuousConfigDoFValues > | m_global_dofs |
std::map< std::string, LocalContinuousConfigDoFValues > | m_local_dofs |
double | m_tol |
CASM::ConfigDoF::ConfigDoF | ( | Index | _N_sublat, |
Index | _N_vol, | ||
GlobalInfoContainerType const & | global_dof_info, | ||
LocalInfoContainerType const & | local_dof_info, | ||
std::vector< SymGroupRepID > const & | occ_symrep_IDs, | ||
double | _tol | ||
) |
Initialize with number of sites, and dimensionality of global and local DoFs GlobalInfoContainerType is an iterable container having value_type std::pair<DoFKey,ContinuousDoFInfo> LocalInfoContainerType is an iterable container having value_type std::pair<DoFKey,std::vector<ContinuousDoFInfo> > OccInfoContainerType is an iterable container having value_type SymGroupRepID
ConfigDoF constructor
_N_sublat | Number of sublattices in corresponding prim |
_N_vol | Supercell volume, as multiples of the prim volume |
global_dof_info | GlobalInfoContainerType is an iterable container of value_type std::pair<DoFKey,ContinuousDoFInfo>. |
local_dof_info | LocalInfoContainerType is an iterable container of value_type std::pair<DoFKey,std::vector<ContinuousDoFInfo> > Note: Typically Structure has already been included when a ConfigDoF is constructed, in which case it is best to use make_configdof from clex/ConfigDoFTools.hh. |
Definition at line 24 of file ConfigDoF_impl.hh.
ConfigDoF & CASM::ConfigDoF::apply_sym | ( | PermuteIterator const & | it | ) |
Update DoF values using the effect of symmetry, including permutation among sites
Update DoF values using the effect of symmetry, including permutation among sites
The effect is: *this = permute_iterator * (*this)
Definition at line 182 of file ConfigDoF.cc.
Update DoF values using only the effect of symmetry on the value at each site, without permutation among sites
Definition at line 220 of file ConfigDoF.cc.
GlobalContinuousConfigDoFValues & CASM::ConfigDoF::global_dof | ( | DoFKey const & | _key | ) |
Definition at line 79 of file ConfigDoF.cc.
GlobalContinuousConfigDoFValues const & CASM::ConfigDoF::global_dof | ( | DoFKey const & | _key | ) | const |
Definition at line 69 of file ConfigDoF.cc.
std::map< DoFKey, GlobalContinuousConfigDoFValues > const & CASM::ConfigDoF::global_dofs | ( | ) | const |
Definition at line 65 of file ConfigDoF.cc.
bool CASM::ConfigDoF::has_global_dof | ( | DoFKey const & | _key | ) | const |
Definition at line 88 of file ConfigDoF.cc.
bool CASM::ConfigDoF::has_local_dof | ( | DoFKey const & | _key | ) | const |
Definition at line 143 of file ConfigDoF.cc.
bool CASM::ConfigDoF::has_occupation | ( | ) | const |
Definition at line 60 of file ConfigDoF.cc.
LocalContinuousConfigDoFValues & CASM::ConfigDoF::local_dof | ( | DoFKey const & | _key | ) |
Definition at line 134 of file ConfigDoF.cc.
LocalContinuousConfigDoFValues const & CASM::ConfigDoF::local_dof | ( | DoFKey const & | _key | ) | const |
Definition at line 124 of file ConfigDoF.cc.
std::map< DoFKey, LocalContinuousConfigDoFValues > const & CASM::ConfigDoF::local_dofs | ( | ) | const |
Definition at line 119 of file ConfigDoF.cc.
Index CASM::ConfigDoF::n_sublat | ( | ) | const |
Number of sublattices in ConfigDoF.
Definition at line 17 of file ConfigDoF.cc.
Index CASM::ConfigDoF::n_vol | ( | ) | const |
Integer volume of ConfigDoF.
Definition at line 14 of file ConfigDoF.cc.
int & CASM::ConfigDoF::occ | ( | Index | i | ) |
Reference occupation value on site i.
Definition at line 34 of file ConfigDoF.cc.
const int & CASM::ConfigDoF::occ | ( | Index | i | ) | const |
Const reference to occupation value on site i.
Definition at line 37 of file ConfigDoF.cc.
Eigen::VectorXi const & CASM::ConfigDoF::occupation | ( | ) | const |
Const reference occupation values.
Definition at line 56 of file ConfigDoF.cc.
void CASM::ConfigDoF::set_global_dof | ( | DoFKey const & | _key, |
Eigen::Ref< const Eigen::VectorXd > const & | _val | ||
) |
Set global continuous DoF values.
Set global continuous DoF values
std::runtime_error | ("Attempting to assign global ConfigDoF values...") if "_key" does not exist in global dofs |
std::runtime_error | ("Size mismatch in ConfigDoF::set_global_dof...") if input vector size does not match current size |
Definition at line 99 of file ConfigDoF.cc.
void CASM::ConfigDoF::set_local_dof | ( | DoFKey const & | _key, |
Eigen::Ref< const Eigen::MatrixXd > const & | _val | ||
) |
Set local continuous DoF values.
Set local continuous DoF values
std::runtime_error | ("Attempting to assign local ConfigDoF values...") if "_key" does not exist in local dofs |
std::runtime_error | ("Size mismatch in ConfigDoF::set_local_dof...") if input matrix size does not match current matrix size |
Definition at line 154 of file ConfigDoF.cc.
void CASM::ConfigDoF::set_occupation | ( | Eigen::Ref< const Eigen::VectorXi > const & | _occupation | ) |
Set occupation values.
Set occupation values
std::runtime_error | ("Size mismatch in ConfigDoF::set_occupation...") if input vector size does not match current size |
Definition at line 43 of file ConfigDoF.cc.
void CASM::ConfigDoF::setZero | ( | ) |
Set all DoF values to zero.
Definition at line 23 of file ConfigDoF.cc.
Index CASM::ConfigDoF::size | ( | ) | const |
Number of sites in the ConfigDoF.
Definition at line 11 of file ConfigDoF.cc.
void CASM::ConfigDoF::swap | ( | ConfigDoF & | RHS | ) |
Definition at line 253 of file ConfigDoF.cc.
double CASM::ConfigDoF::tol | ( | ) | const |
Tolerance for comparison of continuous DoF values.
Definition at line 20 of file ConfigDoF.cc.
|
private |
Definition at line 203 of file ConfigDoF.hh.
|
private |
Definition at line 205 of file ConfigDoF.hh.
|
private |
Definition at line 201 of file ConfigDoF.hh.
|
mutableprivate |
Tolerance used for transformation to canonical form – used also for continuous DoF comparisons, since comparisons are only meaningful to within the tolerance used for finding the canonical form.
Definition at line 210 of file ConfigDoF.hh.