CASM  1.1.0
A Clusters Approach to Statistical Mechanics
StrainDoFTraits.cc
Go to the documentation of this file.
2 
8 #include "casm/clex/ClexBasis.hh"
9 #include "casm/clex/ConfigDoF.hh"
20 
21 namespace CASM {
22 namespace DoF_impl {
23 
28 std::pair<Eigen::MatrixXd, std::set<std::string> > StrainDoFTraits::find_values(
29  std::map<std::string, Eigen::MatrixXd> const &values) const {
30  std::pair<Eigen::MatrixXd, std::set<std::string> > result;
31  for (auto const &val : values) {
32  std::string valname = AnisoValTraits::name_suffix(val.first);
33  auto pos = valname.find("strain");
34  if (pos != std::string::npos) {
35  StrainConverter c(valname.substr(0, pos));
37  c.set_mode(m_metric);
38 
39  result.first = c.unrolled_strain_metric(F);
40  result.second.insert(val.first);
41  return result;
42  }
43  }
44  throw std::runtime_error(
45  "Could not identify DoF values for DoF '" + (this->name()) +
46  "' from provided list of tabulated structure properties.");
47  return result;
48 }
49 
53  Structure const &_prim,
54  std::vector<Orbit<PrimPeriodicSymCompare<IntegralCluster> > > &_asym_unit,
55  BasisFunctionSpecs const &_basis_function_specs) const {
56  // std::cout << "Using " << func_type << " site basis functions." << std::endl
57  // << std::endl;
58  if (_prim.structure().global_dofs().find(name()) ==
59  _prim.structure().global_dofs().end())
60  return std::vector<BasisSet>();
61 
62  std::vector<BasisSet> result(1);
63  result[0].set_variable_basis(adapter::Adapter<CASM::DoFSet, xtal::DoFSet>()(
64  _prim.structure().global_dof(name()),
65  _prim.global_dof_symrep_ID(name())));
66 
67  return result;
68 }
69 
72  BasicStructure const &_reference,
73  SimpleStructure &_struc) const {
74  Eigen::VectorXd unrolled_metric = _dof.global_dof(name()).standard_values();
76  Eigen::Matrix3d F = c.unrolled_strain_metric_to_F(unrolled_metric);
77 
78  _struc.lat_column_mat = F * _struc.lat_column_mat;
79  _struc.mol_info.coords = F * _struc.mol_info.coords;
80  _struc.atom_info.coords = F * _struc.atom_info.coords;
81 }
82 
84  ConfigDoF const &_dof, BasicStructure const &_reference) const {
85  Eigen::VectorXd unrolled_metric = _dof.global_dof(name()).standard_values();
87  Eigen::Matrix3d F = c.unrolled_strain_metric_to_F(unrolled_metric);
88 
89  jsonParser json;
90  json["deformation"] = F;
91  to_json_array(unrolled_metric, json["value"]);
92  return json;
93 }
94 
95 } // namespace DoF_impl
96 
97 namespace DoFType {
99 
101 
103 
104 } // namespace DoFType
105 } // namespace CASM
static std::string name_suffix(std::string const &_name, char delim='_')
GlobalContinuousConfigDoFValues const & global_dof(DoFKey const &_key) const
Definition: ConfigDoF.cc:69
jsonParser dof_to_json(ConfigDoF const &_dof, BasicStructure const &_reference) const override
Serialize strain DoF values from ConfigDoF.
std::vector< BasisSet > construct_site_bases(Structure const &_prim, std::vector< Orbit< PrimPeriodicSymCompare< IntegralCluster > > > &_asym_unit, BasisFunctionSpecs const &_basis_function_specs) const override
Construct the site basis (if DOF_MODE is LOCAL) for a DoF, given its site.
void apply_dof(ConfigDoF const &_dof, BasicStructure const &_reference, SimpleStructure &_struc) const override
Transforms SimpleSructure.
std::pair< Eigen::MatrixXd, std::set< std::string > > find_values(std::map< std::string, Eigen::MatrixXd > const &values) const override
Retrieve the standard values for a DoF from dictionary of properties from properties....
std::string const & name() const
Definition: DoFTraits.hh:69
Eigen::MatrixXd standard_values() const
Get global DoF values as standard DoF values.
An Orbit of Element.
Definition: Orbit.hh:43
Matrix3d unrolled_strain_metric_to_F(Eigen::Ref< const VectorXd > const &E) const
VectorXd unrolled_strain_metric(Eigen::Ref< const Matrix3d > const &F) const
void set_mode(const std::string &mode_name)
Structure specifies the lattice and atomic basis of a crystal.
Definition: Structure.hh:30
const xtal::BasicStructure & structure() const
Definition: Structure.hh:92
SymGroupRepID global_dof_symrep_ID(const std::string dof_name) const
Definition: Structure.cc:137
BasicStructure specifies the lattice and atomic basis of a crystal.
std::map< DoFKey, DoFSet > const & global_dofs() const
DoFSet const & global_dof(std::string const &dof_type) const
Representation of a crystal of molecular and/or atomic occupants, and any additional properties....
CASM::jsonParser & to_json_array(const Eigen::MatrixBase< Derived > &value, CASM::jsonParser &json)
Write Eigen Matrix with 1 row or 1 column to JSON array.
Definition: json_io.hh:313
DoF_impl::StrainDoFTraits GLstrain()
DoF_impl::StrainDoFTraits Hstrain()
DoF_impl::StrainDoFTraits EAstrain()
Main CASM namespace.
Definition: APICommand.hh:8
Eigen::Matrix3d Matrix3d
Eigen::VectorXd VectorXd
Specify how to construct basis functions.
Eigen::MatrixXd coords
(3 x names.size()) matrix of coordinates. coords.col(i) is Cartesian coordinate of site 'i'